Remove constant acceleration groups
[alexxy/gromacs.git] / src / gromacs / mdlib / update.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "update.h"
41
42 #include <cmath>
43 #include <cstdio>
44
45 #include <algorithm>
46 #include <memory>
47
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/gmxlib/nrnb.h"
52 #include "gromacs/listed_forces/disre.h"
53 #include "gromacs/listed_forces/orires.h"
54 #include "gromacs/math/functions.h"
55 #include "gromacs/math/invertmatrix.h"
56 #include "gromacs/math/paddedvector.h"
57 #include "gromacs/math/units.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/math/vecdump.h"
60 #include "gromacs/mdlib/boxdeformation.h"
61 #include "gromacs/mdlib/constr.h"
62 #include "gromacs/mdlib/gmx_omp_nthreads.h"
63 #include "gromacs/mdlib/mdatoms.h"
64 #include "gromacs/mdlib/stat.h"
65 #include "gromacs/mdlib/tgroup.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/fcdata.h"
68 #include "gromacs/mdtypes/group.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/mdtypes/mdatom.h"
72 #include "gromacs/mdtypes/state.h"
73 #include "gromacs/pbcutil/boxutilities.h"
74 #include "gromacs/pbcutil/pbc.h"
75 #include "gromacs/pulling/pull.h"
76 #include "gromacs/random/tabulatednormaldistribution.h"
77 #include "gromacs/random/threefry.h"
78 #include "gromacs/simd/simd.h"
79 #include "gromacs/timing/wallcycle.h"
80 #include "gromacs/topology/atoms.h"
81 #include "gromacs/utility/exceptions.h"
82 #include "gromacs/utility/fatalerror.h"
83 #include "gromacs/utility/futil.h"
84 #include "gromacs/utility/gmxassert.h"
85 #include "gromacs/utility/gmxomp.h"
86 #include "gromacs/utility/smalloc.h"
87
88 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
89
90 struct gmx_sd_const_t
91 {
92     double em = 0;
93 };
94
95 struct gmx_sd_sigma_t
96 {
97     real V = 0;
98 };
99
100 struct gmx_stochd_t
101 {
102     /* BD stuff */
103     std::vector<real> bd_rf;
104     /* SD stuff */
105     std::vector<gmx_sd_const_t> sdc;
106     std::vector<gmx_sd_sigma_t> sdsig;
107     /* andersen temperature control stuff */
108     std::vector<bool> randomize_group;
109     std::vector<real> boltzfac;
110
111     explicit gmx_stochd_t(const t_inputrec& inputRecord);
112 };
113
114 //! pImpled implementation for Update
115 class Update::Impl
116 {
117 public:
118     //! Constructor
119     Impl(const t_inputrec& inputRecord, BoxDeformation* boxDeformation);
120     //! Destructor
121     ~Impl() = default;
122
123     void update_coords(const t_inputrec&                                inputRecord,
124                        int64_t                                          step,
125                        const t_mdatoms*                                 md,
126                        t_state*                                         state,
127                        const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
128                        const t_fcdata&                                  fcdata,
129                        const gmx_ekindata_t*                            ekind,
130                        const matrix                                     M,
131                        int                                              UpdatePart,
132                        const t_commrec*                                 cr,
133                        bool                                             haveConstraints);
134
135     void finish_update(const t_inputrec& inputRecord,
136                        const t_mdatoms*  md,
137                        t_state*          state,
138                        gmx_wallcycle_t   wcycle,
139                        bool              haveConstraints);
140
141     void update_sd_second_half(const t_inputrec& inputRecord,
142                                int64_t           step,
143                                real*             dvdlambda,
144                                const t_mdatoms*  md,
145                                t_state*          state,
146                                const t_commrec*  cr,
147                                t_nrnb*           nrnb,
148                                gmx_wallcycle_t   wcycle,
149                                gmx::Constraints* constr,
150                                bool              do_log,
151                                bool              do_ene);
152
153     void update_for_constraint_virial(const t_inputrec&                                inputRecord,
154                                       const t_mdatoms&                                 md,
155                                       const t_state&                                   state,
156                                       const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
157                                       const gmx_ekindata_t&                            ekind);
158
159     void update_temperature_constants(const t_inputrec& inputRecord);
160
161     const std::vector<bool>& getAndersenRandomizeGroup() const { return sd_.randomize_group; }
162
163     const std::vector<real>& getBoltzmanFactor() const { return sd_.boltzfac; }
164
165     PaddedVector<RVec>* xp() { return &xp_; }
166
167     BoxDeformation* deform() const { return deform_; }
168
169 private:
170     //! stochastic dynamics struct
171     gmx_stochd_t sd_;
172     //! xprime for constraint algorithms
173     PaddedVector<RVec> xp_;
174     //! Box deformation handler (or nullptr if inactive).
175     BoxDeformation* deform_ = nullptr;
176 };
177
178 Update::Update(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
179     impl_(new Impl(inputRecord, boxDeformation)){};
180
181 Update::~Update(){};
182
183 const std::vector<bool>& Update::getAndersenRandomizeGroup() const
184 {
185     return impl_->getAndersenRandomizeGroup();
186 }
187
188 const std::vector<real>& Update::getBoltzmanFactor() const
189 {
190     return impl_->getBoltzmanFactor();
191 }
192
193 PaddedVector<RVec>* Update::xp()
194 {
195     return impl_->xp();
196 }
197
198 BoxDeformation* Update::deform() const
199 {
200     return impl_->deform();
201 }
202
203 void Update::update_coords(const t_inputrec&                                inputRecord,
204                            int64_t                                          step,
205                            const t_mdatoms*                                 md,
206                            t_state*                                         state,
207                            const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
208                            const t_fcdata&                                  fcdata,
209                            const gmx_ekindata_t*                            ekind,
210                            const matrix                                     M,
211                            int                                              updatePart,
212                            const t_commrec*                                 cr,
213                            const bool                                       haveConstraints)
214 {
215     return impl_->update_coords(
216             inputRecord, step, md, state, f, fcdata, ekind, M, updatePart, cr, haveConstraints);
217 }
218
219 void Update::finish_update(const t_inputrec& inputRecord,
220                            const t_mdatoms*  md,
221                            t_state*          state,
222                            gmx_wallcycle_t   wcycle,
223                            const bool        haveConstraints)
224 {
225     return impl_->finish_update(inputRecord, md, state, wcycle, haveConstraints);
226 }
227
228 void Update::update_sd_second_half(const t_inputrec& inputRecord,
229                                    int64_t           step,
230                                    real*             dvdlambda,
231                                    const t_mdatoms*  md,
232                                    t_state*          state,
233                                    const t_commrec*  cr,
234                                    t_nrnb*           nrnb,
235                                    gmx_wallcycle_t   wcycle,
236                                    gmx::Constraints* constr,
237                                    bool              do_log,
238                                    bool              do_ene)
239 {
240     return impl_->update_sd_second_half(
241             inputRecord, step, dvdlambda, md, state, cr, nrnb, wcycle, constr, do_log, do_ene);
242 }
243
244 void Update::update_for_constraint_virial(const t_inputrec& inputRecord,
245                                           const t_mdatoms&  md,
246                                           const t_state&    state,
247                                           const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
248                                           const gmx_ekindata_t&                            ekind)
249 {
250     return impl_->update_for_constraint_virial(inputRecord, md, state, f, ekind);
251 }
252
253 void Update::update_temperature_constants(const t_inputrec& inputRecord)
254 {
255     return impl_->update_temperature_constants(inputRecord);
256 }
257
258 /*! \brief Sets the velocities of virtual sites to zero */
259 static void clearVsiteVelocities(int start, int nrend, const unsigned short* particleType, rvec* gmx_restrict v)
260 {
261     for (int a = start; a < nrend; a++)
262     {
263         if (particleType[a] == eptVSite)
264         {
265             clear_rvec(v[a]);
266         }
267     }
268 }
269
270 /*! \brief Sets whether we store the updated velocities */
271 enum class StoreUpdatedVelocities
272 {
273     yes, //!< Store the updated velocities
274     no   //!< Do not store the updated velocities
275 };
276
277 /*! \brief Sets the number of different temperature coupling values */
278 enum class NumTempScaleValues
279 {
280     single,  //!< Single T-scaling value (either one group or all values =1)
281     multiple //!< Multiple T-scaling values, need to use T-group indices
282 };
283
284 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
285  *
286  * Note that this enum is only used in updateMDLeapfrogSimple(), which does
287  * not handle fully anistropic Parrinello-Rahman scaling, so we only have
288  * options \p no and \p diagonal here and no anistropic option.
289  */
290 enum class ApplyParrinelloRahmanVScaling
291 {
292     no,      //!< Do not apply velocity scaling (not a PR-coupling run or step)
293     diagonal //!< Apply velocity scaling using a diagonal matrix
294 };
295
296 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
297  *
298  * \tparam       storeUpdatedVelocities Tells whether we should store the updated velocities
299  * \tparam       numTempScaleValues     The number of different T-couple values
300  * \tparam       applyPRVScaling        Apply Parrinello-Rahman velocity scaling
301  * \param[in]    start                  Index of first atom to update
302  * \param[in]    nrend                  Last atom to update: \p nrend - 1
303  * \param[in]    dt                     The time step
304  * \param[in]    dtPressureCouple       Time step for pressure coupling
305  * \param[in]    invMassPerDim          1/mass per atom and dimension
306  * \param[in]    tcstat                 Temperature coupling information
307  * \param[in]    cTC                    T-coupling group index per atom
308  * \param[in]    pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
309  * \param[in]    x                      Input coordinates
310  * \param[out]   xprime                 Updated coordinates
311  * \param[inout] v                      Velocities
312  * \param[in]    f                      Forces
313  *
314  * We expect this template to get good SIMD acceleration by most compilers,
315  * unlike the more complex general template.
316  * Note that we might get even better SIMD acceleration when we introduce
317  * aligned (and padded) memory, possibly with some hints for the compilers.
318  */
319 template<StoreUpdatedVelocities storeUpdatedVelocities, NumTempScaleValues numTempScaleValues, ApplyParrinelloRahmanVScaling applyPRVScaling>
320 static void updateMDLeapfrogSimple(int         start,
321                                    int         nrend,
322                                    real        dt,
323                                    real        dtPressureCouple,
324                                    const rvec* gmx_restrict          invMassPerDim,
325                                    gmx::ArrayRef<const t_grp_tcstat> tcstat,
326                                    const unsigned short*             cTC,
327                                    const rvec                        pRVScaleMatrixDiagonal,
328                                    const rvec* gmx_restrict x,
329                                    rvec* gmx_restrict xprime,
330                                    rvec* gmx_restrict v,
331                                    const rvec* gmx_restrict f)
332 {
333     real lambdaGroup;
334
335     if (numTempScaleValues == NumTempScaleValues::single)
336     {
337         lambdaGroup = tcstat[0].lambda;
338     }
339
340     for (int a = start; a < nrend; a++)
341     {
342         if (numTempScaleValues == NumTempScaleValues::multiple)
343         {
344             lambdaGroup = tcstat[cTC[a]].lambda;
345         }
346
347         for (int d = 0; d < DIM; d++)
348         {
349             /* Note that using rvec invMassPerDim results in more efficient
350              * SIMD code, but this increases the cache pressure.
351              * For large systems with PME on the CPU this slows down the
352              * (then already slow) update by 20%. If all data remains in cache,
353              * using rvec is much faster.
354              */
355             real vNew = lambdaGroup * v[a][d] + f[a][d] * invMassPerDim[a][d] * dt;
356
357             if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
358             {
359                 vNew -= dtPressureCouple * pRVScaleMatrixDiagonal[d] * v[a][d];
360             }
361             if (storeUpdatedVelocities == StoreUpdatedVelocities::yes)
362             {
363                 v[a][d] = vNew;
364             }
365             xprime[a][d] = x[a][d] + vNew * dt;
366         }
367     }
368 }
369
370 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
371 #    define GMX_HAVE_SIMD_UPDATE 1
372 #else
373 #    define GMX_HAVE_SIMD_UPDATE 0
374 #endif
375
376 #if GMX_HAVE_SIMD_UPDATE
377
378 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
379  *
380  * The loaded output is:
381  * \p r0: { r[index][XX], r[index][YY], ... }
382  * \p r1: { ... }
383  * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
384  *
385  * \param[in]  r      Real to an rvec array, has to be aligned to SIMD register width
386  * \param[in]  index  Index of the first rvec triplet of reals to load
387  * \param[out] r0     Pointer to first SIMD register
388  * \param[out] r1     Pointer to second SIMD register
389  * \param[out] r2     Pointer to third SIMD register
390  */
391 static inline void simdLoadRvecs(const rvec* r, int index, SimdReal* r0, SimdReal* r1, SimdReal* r2)
392 {
393     const real* realPtr = r[index];
394
395     GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
396
397     *r0 = simdLoad(realPtr + 0 * GMX_SIMD_REAL_WIDTH);
398     *r1 = simdLoad(realPtr + 1 * GMX_SIMD_REAL_WIDTH);
399     *r2 = simdLoad(realPtr + 2 * GMX_SIMD_REAL_WIDTH);
400 }
401
402 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
403  *
404  * The stored output is:
405  * \p r[index] = { { r0[0], r0[1], ... }
406  * ...
407  * \p r[index+GMX_SIMD_REAL_WIDTH-1] =  { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
408  *
409  * \param[out] r      Pointer to an rvec array, has to be aligned to SIMD register width
410  * \param[in]  index  Index of the first rvec triplet of reals to store to
411  * \param[in]  r0     First SIMD register
412  * \param[in]  r1     Second SIMD register
413  * \param[in]  r2     Third SIMD register
414  */
415 static inline void simdStoreRvecs(rvec* r, int index, SimdReal r0, SimdReal r1, SimdReal r2)
416 {
417     real* realPtr = r[index];
418
419     GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
420
421     store(realPtr + 0 * GMX_SIMD_REAL_WIDTH, r0);
422     store(realPtr + 1 * GMX_SIMD_REAL_WIDTH, r1);
423     store(realPtr + 2 * GMX_SIMD_REAL_WIDTH, r2);
424 }
425
426 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
427  *
428  * \tparam       storeUpdatedVelocities Tells whether we should store the updated velocities
429  * \param[in]    start                  Index of first atom to update
430  * \param[in]    nrend                  Last atom to update: \p nrend - 1
431  * \param[in]    dt                     The time step
432  * \param[in]    invMass                1/mass per atom
433  * \param[in]    tcstat                 Temperature coupling information
434  * \param[in]    x                      Input coordinates
435  * \param[out]   xprime                 Updated coordinates
436  * \param[inout] v                      Velocities
437  * \param[in]    f                      Forces
438  */
439 template<StoreUpdatedVelocities storeUpdatedVelocities>
440 static void updateMDLeapfrogSimpleSimd(int         start,
441                                        int         nrend,
442                                        real        dt,
443                                        const real* gmx_restrict          invMass,
444                                        gmx::ArrayRef<const t_grp_tcstat> tcstat,
445                                        const rvec* gmx_restrict x,
446                                        rvec* gmx_restrict xprime,
447                                        rvec* gmx_restrict v,
448                                        const rvec* gmx_restrict f)
449 {
450     SimdReal timestep(dt);
451     SimdReal lambdaSystem(tcstat[0].lambda);
452
453     /* We declare variables here, since code is often slower when declaring them inside the loop */
454
455     /* Note: We should implement a proper PaddedVector, so we don't need this check */
456     GMX_ASSERT(isSimdAligned(invMass), "invMass should be aligned");
457
458     for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
459     {
460         SimdReal invMass0, invMass1, invMass2;
461         expandScalarsToTriplets(simdLoad(invMass + a), &invMass0, &invMass1, &invMass2);
462
463         SimdReal v0, v1, v2;
464         SimdReal f0, f1, f2;
465         simdLoadRvecs(v, a, &v0, &v1, &v2);
466         simdLoadRvecs(f, a, &f0, &f1, &f2);
467
468         v0 = fma(f0 * invMass0, timestep, lambdaSystem * v0);
469         v1 = fma(f1 * invMass1, timestep, lambdaSystem * v1);
470         v2 = fma(f2 * invMass2, timestep, lambdaSystem * v2);
471
472         if (storeUpdatedVelocities == StoreUpdatedVelocities::yes)
473         {
474             simdStoreRvecs(v, a, v0, v1, v2);
475         }
476
477         SimdReal x0, x1, x2;
478         simdLoadRvecs(x, a, &x0, &x1, &x2);
479
480         SimdReal xprime0 = fma(v0, timestep, x0);
481         SimdReal xprime1 = fma(v1, timestep, x1);
482         SimdReal xprime2 = fma(v2, timestep, x2);
483
484         simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
485     }
486 }
487
488 #endif // GMX_HAVE_SIMD_UPDATE
489
490 /*! \brief Sets the NEMD acceleration type */
491 enum class AccelerationType
492 {
493     none,
494     cosine
495 };
496
497 /*! \brief Integrate using leap-frog with support for everything.
498  *
499  * \tparam        accelerationType  Type of NEMD acceleration.
500  * \param[in]     start             Index of first atom to update.
501  * \param[in]     nrend             Last atom to update: \p nrend - 1.
502  * \param[in]     doNoseHoover      If to apply Nose-Hoover T-coupling.
503  * \param[in]     dt                The time step.
504  * \param[in]     dtPressureCouple  Time step for pressure coupling, is 0 when no pressure
505  *                                  coupling should be applied at this step.
506  * \param[in]     md                Atom properties.
507  * \param[in]     ekind             Kinetic energy data.
508  * \param[in]     box               The box dimensions.
509  * \param[in]     x                 Input coordinates.
510  * \param[out]    xprime            Updated coordinates.
511  * \param[inout]  v                 Velocities.
512  * \param[in]     f                 Forces.
513  * \param[in]     nh_vxi            Nose-Hoover velocity scaling factors.
514  * \param[in]     nsttcouple        Frequency of the temperature coupling steps.
515  * \param[in]     M                 Parrinello-Rahman scaling matrix.
516  */
517 template<AccelerationType accelerationType>
518 static void updateMDLeapfrogGeneral(int                   start,
519                                     int                   nrend,
520                                     bool                  doNoseHoover,
521                                     real                  dt,
522                                     real                  dtPressureCouple,
523                                     const t_mdatoms*      md,
524                                     const gmx_ekindata_t* ekind,
525                                     const matrix          box,
526                                     const rvec* gmx_restrict x,
527                                     rvec* gmx_restrict xprime,
528                                     rvec* gmx_restrict v,
529                                     const rvec* gmx_restrict f,
530                                     const double* gmx_restrict nh_vxi,
531                                     const int                  nsttcouple,
532                                     const matrix               M)
533 {
534     /* This is a version of the leap-frog integrator that supports
535      * all combinations of T-coupling, P-coupling and NEMD.
536      * Nose-Hoover uses the reversible leap-frog integrator from
537      * Holian et al. Phys Rev E 52(3) : 2338, 1995
538      */
539
540     gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
541     const unsigned short*             cTC    = md->cTC;
542
543     const rvec* gmx_restrict invMassPerDim = md->invMassPerDim;
544
545     /* Initialize group values, changed later when multiple groups are used */
546     int gt = 0;
547
548     real omega_Z = 2 * static_cast<real>(M_PI) / box[ZZ][ZZ];
549
550     for (int n = start; n < nrend; n++)
551     {
552         if (cTC)
553         {
554             gt = cTC[n];
555         }
556         real lg = tcstat[gt].lambda;
557
558         rvec vRel;
559         real cosineZ, vCosine;
560 #ifdef __INTEL_COMPILER
561 #    pragma warning(disable : 280)
562 #endif
563         switch (accelerationType)
564         {
565             case AccelerationType::none: copy_rvec(v[n], vRel); break;
566             case AccelerationType::cosine:
567                 cosineZ = std::cos(x[n][ZZ] * omega_Z);
568                 vCosine = cosineZ * ekind->cosacc.vcos;
569                 /* Avoid scaling the cosine profile velocity */
570                 copy_rvec(v[n], vRel);
571                 vRel[XX] -= vCosine;
572                 break;
573         }
574
575         real factorNH = 0.0;
576         if (doNoseHoover)
577         {
578             /* Here we account for multiple time stepping, by increasing
579              * the Nose-Hoover correction by nsttcouple
580              * TODO: This can be pre-computed.
581              */
582             factorNH = 0.5 * nsttcouple * dt * nh_vxi[gt];
583         }
584
585         for (int d = 0; d < DIM; d++)
586         {
587             real vNew = (lg * vRel[d]
588                          + (f[n][d] * invMassPerDim[n][d] * dt - factorNH * vRel[d]
589                             - dtPressureCouple * iprod(M[d], vRel)))
590                         / (1 + factorNH);
591             switch (accelerationType)
592             {
593                 case AccelerationType::none: break;
594                 case AccelerationType::cosine:
595                     if (d == XX)
596                     {
597                         /* Add back the mean velocity and apply acceleration */
598                         vNew += vCosine + cosineZ * ekind->cosacc.cos_accel * dt;
599                     }
600                     break;
601             }
602             v[n][d]      = vNew;
603             xprime[n][d] = x[n][d] + vNew * dt;
604         }
605     }
606 }
607
608 /*! \brief Handles the Leap-frog MD x and v integration */
609 static void do_update_md(int         start,
610                          int         nrend,
611                          real        dt,
612                          int64_t     step,
613                          const rvec* gmx_restrict x,
614                          rvec* gmx_restrict xprime,
615                          rvec* gmx_restrict v,
616                          const rvec* gmx_restrict f,
617                          const int                etc,
618                          const int                epc,
619                          const int                nsttcouple,
620                          const int                nstpcouple,
621                          const t_mdatoms*         md,
622                          const gmx_ekindata_t*    ekind,
623                          const matrix             box,
624                          const double* gmx_restrict nh_vxi,
625                          const matrix               M)
626 {
627     GMX_ASSERT(nrend == start || xprime != x,
628                "For SIMD optimization certain compilers need to have xprime != x");
629
630     /* Note: Berendsen pressure scaling is handled after do_update_md() */
631     bool doTempCouple = (etc != etcNO && do_per_step(step + nsttcouple - 1, nsttcouple));
632     bool doNoseHoover = (etc == etcNOSEHOOVER && doTempCouple);
633     bool doParrinelloRahman =
634             (epc == epcPARRINELLORAHMAN && do_per_step(step + nstpcouple - 1, nstpcouple));
635     bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
636
637     real dtPressureCouple = (doParrinelloRahman ? nstpcouple * dt : 0);
638
639     /* NEMD cosine acceleration is applied in updateMDLeapFrogGeneral */
640     bool doAcceleration = (ekind->cosacc.cos_accel != 0);
641
642     if (doNoseHoover || doPROffDiagonal || doAcceleration)
643     {
644         matrix stepM;
645         if (!doParrinelloRahman)
646         {
647             /* We should not apply PR scaling at this step */
648             clear_mat(stepM);
649         }
650         else
651         {
652             copy_mat(M, stepM);
653         }
654
655         if (!doAcceleration)
656         {
657             updateMDLeapfrogGeneral<AccelerationType::none>(
658                     start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
659         }
660         else
661         {
662             updateMDLeapfrogGeneral<AccelerationType::cosine>(
663                     start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
664         }
665     }
666     else
667     {
668         /* Use a simple and thus more efficient integration loop. */
669         /* The simple loop does not check for particle type (so it can
670          * be vectorized), which means we need to clear the velocities
671          * of virtual sites in advance, when present. Note that vsite
672          * velocities are computed after update and constraints from
673          * their displacement.
674          */
675         if (md->haveVsites)
676         {
677             /* Note: The overhead of this loop is completely neligible */
678             clearVsiteVelocities(start, nrend, md->ptype, v);
679         }
680
681         /* We determine if we have a single T-coupling lambda value for all
682          * atoms. That allows for better SIMD acceleration in the template.
683          * If we do not do temperature coupling (in the run or this step),
684          * all scaling values are 1, so we effectively have a single value.
685          */
686         bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
687
688         /* Extract some pointers needed by all cases */
689         const unsigned short*             cTC           = md->cTC;
690         gmx::ArrayRef<const t_grp_tcstat> tcstat        = ekind->tcstat;
691         const rvec*                       invMassPerDim = md->invMassPerDim;
692
693         if (doParrinelloRahman)
694         {
695             GMX_ASSERT(!doPROffDiagonal,
696                        "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling "
697                        "matrices");
698
699             rvec diagM;
700             for (int d = 0; d < DIM; d++)
701             {
702                 diagM[d] = M[d][d];
703             }
704
705             if (haveSingleTempScaleValue)
706             {
707                 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::diagonal>(
708                         start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
709             }
710             else
711             {
712                 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::diagonal>(
713                         start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
714             }
715         }
716         else
717         {
718             if (haveSingleTempScaleValue)
719             {
720                 /* Note that modern compilers are pretty good at vectorizing
721                  * updateMDLeapfrogSimple(). But the SIMD version will still
722                  * be faster because invMass lowers the cache pressure
723                  * compared to invMassPerDim.
724                  */
725 #if GMX_HAVE_SIMD_UPDATE
726                 /* Check if we can use invmass instead of invMassPerDim */
727                 if (!md->havePartiallyFrozenAtoms)
728                 {
729                     updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::yes>(
730                             start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
731                 }
732                 else
733 #endif
734                 {
735                     updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
736                             start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
737                 }
738             }
739             else
740             {
741                 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::no>(
742                         start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
743             }
744         }
745     }
746 }
747 /*! \brief Handles the Leap-frog MD x and v integration */
748 static void doUpdateMDDoNotUpdateVelocities(int         start,
749                                             int         nrend,
750                                             real        dt,
751                                             const rvec* gmx_restrict x,
752                                             rvec* gmx_restrict xprime,
753                                             rvec* gmx_restrict v,
754                                             const rvec* gmx_restrict f,
755                                             const t_mdatoms&         md,
756                                             const gmx_ekindata_t&    ekind)
757 {
758     GMX_ASSERT(nrend == start || xprime != x,
759                "For SIMD optimization certain compilers need to have xprime != x");
760
761     gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind.tcstat;
762
763     /* Check if we can use invmass instead of invMassPerDim */
764 #if GMX_HAVE_SIMD_UPDATE
765     if (!md.havePartiallyFrozenAtoms)
766     {
767         updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::no>(
768                 start, nrend, dt, md.invmass, tcstat, x, xprime, v, f);
769     }
770     else
771 #endif
772     {
773         updateMDLeapfrogSimple<StoreUpdatedVelocities::no, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
774                 start, nrend, dt, dt, md.invMassPerDim, tcstat, nullptr, nullptr, x, xprime, v, f);
775     }
776 }
777
778 static void do_update_vv_vel(int                  start,
779                              int                  nrend,
780                              real                 dt,
781                              const ivec           nFreeze[],
782                              const real           invmass[],
783                              const unsigned short ptype[],
784                              const unsigned short cFREEZE[],
785                              rvec                 v[],
786                              const rvec           f[],
787                              gmx_bool             bExtended,
788                              real                 veta,
789                              real                 alpha)
790 {
791     int  gf = 0;
792     int  n, d;
793     real g, mv1, mv2;
794
795     if (bExtended)
796     {
797         g   = 0.25 * dt * veta * alpha;
798         mv1 = std::exp(-g);
799         mv2 = gmx::series_sinhx(g);
800     }
801     else
802     {
803         mv1 = 1.0;
804         mv2 = 1.0;
805     }
806     for (n = start; n < nrend; n++)
807     {
808         real w_dt = invmass[n] * dt;
809         if (cFREEZE)
810         {
811             gf = cFREEZE[n];
812         }
813
814         for (d = 0; d < DIM; d++)
815         {
816             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
817             {
818                 v[n][d] = mv1 * (mv1 * v[n][d] + 0.5 * (w_dt * mv2 * f[n][d]));
819             }
820             else
821             {
822                 v[n][d] = 0.0;
823             }
824         }
825     }
826 } /* do_update_vv_vel */
827
828 static void do_update_vv_pos(int                  start,
829                              int                  nrend,
830                              real                 dt,
831                              const ivec           nFreeze[],
832                              const unsigned short ptype[],
833                              const unsigned short cFREEZE[],
834                              const rvec           x[],
835                              rvec                 xprime[],
836                              const rvec           v[],
837                              gmx_bool             bExtended,
838                              real                 veta)
839 {
840     int  gf = 0;
841     int  n, d;
842     real g, mr1, mr2;
843
844     /* Would it make more sense if Parrinello-Rahman was put here? */
845     if (bExtended)
846     {
847         g   = 0.5 * dt * veta;
848         mr1 = std::exp(g);
849         mr2 = gmx::series_sinhx(g);
850     }
851     else
852     {
853         mr1 = 1.0;
854         mr2 = 1.0;
855     }
856
857     for (n = start; n < nrend; n++)
858     {
859
860         if (cFREEZE)
861         {
862             gf = cFREEZE[n];
863         }
864
865         for (d = 0; d < DIM; d++)
866         {
867             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
868             {
869                 xprime[n][d] = mr1 * (mr1 * x[n][d] + mr2 * dt * v[n][d]);
870             }
871             else
872             {
873                 xprime[n][d] = x[n][d];
874             }
875         }
876     }
877 } /* do_update_vv_pos */
878
879 gmx_stochd_t::gmx_stochd_t(const t_inputrec& inputRecord)
880 {
881     const t_grpopts* opts = &inputRecord.opts;
882     const int        ngtc = opts->ngtc;
883
884     if (inputRecord.eI == eiBD)
885     {
886         bd_rf.resize(ngtc);
887     }
888     else if (EI_SD(inputRecord.eI))
889     {
890         sdc.resize(ngtc);
891         sdsig.resize(ngtc);
892
893         for (int gt = 0; gt < ngtc; gt++)
894         {
895             if (opts->tau_t[gt] > 0)
896             {
897                 sdc[gt].em = std::exp(-inputRecord.delta_t / opts->tau_t[gt]);
898             }
899             else
900             {
901                 /* No friction and noise on this group */
902                 sdc[gt].em = 1;
903             }
904         }
905     }
906     else if (ETC_ANDERSEN(inputRecord.etc))
907     {
908         randomize_group.resize(ngtc);
909         boltzfac.resize(ngtc);
910
911         /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
912         /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
913
914         for (int gt = 0; gt < ngtc; gt++)
915         {
916             real reft = std::max<real>(0, opts->ref_t[gt]);
917             if ((opts->tau_t[gt] > 0)
918                 && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
919             {
920                 randomize_group[gt] = true;
921                 boltzfac[gt]        = BOLTZ * opts->ref_t[gt];
922             }
923             else
924             {
925                 randomize_group[gt] = false;
926             }
927         }
928     }
929 }
930
931 void Update::Impl::update_temperature_constants(const t_inputrec& inputRecord)
932 {
933     if (inputRecord.eI == eiBD)
934     {
935         if (inputRecord.bd_fric != 0)
936         {
937             for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
938             {
939                 sd_.bd_rf[gt] = std::sqrt(2.0 * BOLTZ * inputRecord.opts.ref_t[gt]
940                                           / (inputRecord.bd_fric * inputRecord.delta_t));
941             }
942         }
943         else
944         {
945             for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
946             {
947                 sd_.bd_rf[gt] = std::sqrt(2.0 * BOLTZ * inputRecord.opts.ref_t[gt]);
948             }
949         }
950     }
951     if (inputRecord.eI == eiSD1)
952     {
953         for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
954         {
955             real kT = BOLTZ * inputRecord.opts.ref_t[gt];
956             /* The mass is accounted for later, since this differs per atom */
957             sd_.sdsig[gt].V = std::sqrt(kT * (1 - sd_.sdc[gt].em * sd_.sdc[gt].em));
958         }
959     }
960 }
961
962 Update::Impl::Impl(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
963     sd_(inputRecord),
964     deform_(boxDeformation)
965 {
966     update_temperature_constants(inputRecord);
967     xp_.resizeWithPadding(0);
968 }
969
970 void Update::setNumAtoms(int numAtoms)
971 {
972
973     impl_->xp()->resizeWithPadding(numAtoms);
974 }
975
976 /*! \brief Sets the SD update type */
977 enum class SDUpdate : int
978 {
979     ForcesOnly,
980     FrictionAndNoiseOnly,
981     Combined
982 };
983
984 /*! \brief SD integrator update
985  *
986  * Two phases are required in the general case of a constrained
987  * update, the first phase from the contribution of forces, before
988  * applying constraints, and then a second phase applying the friction
989  * and noise, and then further constraining. For details, see
990  * Goga2012.
991  *
992  * Without constraints, the two phases can be combined, for
993  * efficiency.
994  *
995  * Thus three instantiations of this templated function will be made,
996  * two with only one contribution, and one with both contributions. */
997 template<SDUpdate updateType>
998 static void doSDUpdateGeneral(const gmx_stochd_t&  sd,
999                               int                  start,
1000                               int                  nrend,
1001                               real                 dt,
1002                               const ivec           nFreeze[],
1003                               const real           invmass[],
1004                               const unsigned short ptype[],
1005                               const unsigned short cFREEZE[],
1006                               const unsigned short cTC[],
1007                               const rvec           x[],
1008                               rvec                 xprime[],
1009                               rvec                 v[],
1010                               const rvec           f[],
1011                               int64_t              step,
1012                               int                  seed,
1013                               const int*           gatindex)
1014 {
1015     // cTC and cFREEZE can be nullptr any time, but various
1016     // instantiations do not make sense with particular pointer
1017     // values.
1018     if (updateType == SDUpdate::ForcesOnly)
1019     {
1020         GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
1021         GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
1022     }
1023     if (updateType == SDUpdate::FrictionAndNoiseOnly)
1024     {
1025         GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
1026     }
1027     if (updateType == SDUpdate::Combined)
1028     {
1029         GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
1030     }
1031
1032     // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
1033     gmx::ThreeFry2x64<0>                       rng(seed, gmx::RandomDomain::UpdateCoordinates);
1034     gmx::TabulatedNormalDistribution<real, 14> dist;
1035
1036     for (int n = start; n < nrend; n++)
1037     {
1038         int globalAtomIndex = gatindex ? gatindex[n] : n;
1039         rng.restart(step, globalAtomIndex);
1040         dist.reset();
1041
1042         real inverseMass = invmass[n];
1043         real invsqrtMass = std::sqrt(inverseMass);
1044
1045         int freezeGroup      = cFREEZE ? cFREEZE[n] : 0;
1046         int temperatureGroup = cTC ? cTC[n] : 0;
1047
1048         for (int d = 0; d < DIM; d++)
1049         {
1050             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[freezeGroup][d])
1051             {
1052                 if (updateType == SDUpdate::ForcesOnly)
1053                 {
1054                     real vn = v[n][d] + inverseMass * f[n][d] * dt;
1055                     v[n][d] = vn;
1056                     // Simple position update.
1057                     xprime[n][d] = x[n][d] + v[n][d] * dt;
1058                 }
1059                 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
1060                 {
1061                     real vn = v[n][d];
1062                     v[n][d] = (vn * sd.sdc[temperatureGroup].em
1063                                + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1064                     // The previous phase already updated the
1065                     // positions with a full v*dt term that must
1066                     // now be half removed.
1067                     xprime[n][d] = xprime[n][d] + 0.5 * (v[n][d] - vn) * dt;
1068                 }
1069                 else
1070                 {
1071                     real vn = v[n][d] + inverseMass * f[n][d] * dt;
1072                     v[n][d] = (vn * sd.sdc[temperatureGroup].em
1073                                + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1074                     // Here we include half of the friction+noise
1075                     // update of v into the position update.
1076                     xprime[n][d] = x[n][d] + 0.5 * (vn + v[n][d]) * dt;
1077                 }
1078             }
1079             else
1080             {
1081                 // When using constraints, the update is split into
1082                 // two phases, but we only need to zero the update of
1083                 // virtual, shell or frozen particles in at most one
1084                 // of the phases.
1085                 if (updateType != SDUpdate::FrictionAndNoiseOnly)
1086                 {
1087                     v[n][d]      = 0.0;
1088                     xprime[n][d] = x[n][d];
1089                 }
1090             }
1091         }
1092     }
1093 }
1094
1095 static void do_update_sd(int         start,
1096                          int         nrend,
1097                          real        dt,
1098                          int64_t     step,
1099                          const rvec* gmx_restrict x,
1100                          rvec* gmx_restrict xprime,
1101                          rvec* gmx_restrict v,
1102                          const rvec* gmx_restrict f,
1103                          const ivec               nFreeze[],
1104                          const real               invmass[],
1105                          const unsigned short     ptype[],
1106                          const unsigned short     cFREEZE[],
1107                          const unsigned short     cTC[],
1108                          int                      seed,
1109                          const t_commrec*         cr,
1110                          const gmx_stochd_t&      sd,
1111                          bool                     haveConstraints)
1112 {
1113     if (haveConstraints)
1114     {
1115         // With constraints, the SD update is done in 2 parts
1116         doSDUpdateGeneral<SDUpdate::ForcesOnly>(
1117                 sd, start, nrend, dt, nFreeze, invmass, ptype, cFREEZE, nullptr, x, xprime, v, f, step, seed, nullptr);
1118     }
1119     else
1120     {
1121         doSDUpdateGeneral<SDUpdate::Combined>(sd,
1122                                               start,
1123                                               nrend,
1124                                               dt,
1125                                               nFreeze,
1126                                               invmass,
1127                                               ptype,
1128                                               cFREEZE,
1129                                               cTC,
1130                                               x,
1131                                               xprime,
1132                                               v,
1133                                               f,
1134                                               step,
1135                                               seed,
1136                                               DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1137     }
1138 }
1139
1140 static void do_update_bd(int         start,
1141                          int         nrend,
1142                          real        dt,
1143                          int64_t     step,
1144                          const rvec* gmx_restrict x,
1145                          rvec* gmx_restrict xprime,
1146                          rvec* gmx_restrict v,
1147                          const rvec* gmx_restrict f,
1148                          const ivec               nFreeze[],
1149                          const real               invmass[],
1150                          const unsigned short     ptype[],
1151                          const unsigned short     cFREEZE[],
1152                          const unsigned short     cTC[],
1153                          real                     friction_coefficient,
1154                          const real*              rf,
1155                          int                      seed,
1156                          const int*               gatindex)
1157 {
1158     /* note -- these appear to be full step velocities . . .  */
1159     int  gf = 0, gt = 0;
1160     real vn;
1161     real invfr = 0;
1162     int  n, d;
1163     // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1164     // Each 64-bit value is enough for 4 normal distribution table numbers.
1165     gmx::ThreeFry2x64<0>                       rng(seed, gmx::RandomDomain::UpdateCoordinates);
1166     gmx::TabulatedNormalDistribution<real, 14> dist;
1167
1168     if (friction_coefficient != 0)
1169     {
1170         invfr = 1.0 / friction_coefficient;
1171     }
1172
1173     for (n = start; (n < nrend); n++)
1174     {
1175         int ng = gatindex ? gatindex[n] : n;
1176
1177         rng.restart(step, ng);
1178         dist.reset();
1179
1180         if (cFREEZE)
1181         {
1182             gf = cFREEZE[n];
1183         }
1184         if (cTC)
1185         {
1186             gt = cTC[n];
1187         }
1188         for (d = 0; (d < DIM); d++)
1189         {
1190             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
1191             {
1192                 if (friction_coefficient != 0)
1193                 {
1194                     vn = invfr * f[n][d] + rf[gt] * dist(rng);
1195                 }
1196                 else
1197                 {
1198                     /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1199                     vn = 0.5 * invmass[n] * f[n][d] * dt
1200                          + std::sqrt(0.5 * invmass[n]) * rf[gt] * dist(rng);
1201                 }
1202
1203                 v[n][d]      = vn;
1204                 xprime[n][d] = x[n][d] + vn * dt;
1205             }
1206             else
1207             {
1208                 v[n][d]      = 0.0;
1209                 xprime[n][d] = x[n][d];
1210             }
1211         }
1212     }
1213 }
1214
1215 extern void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir)
1216 {
1217     ekinstate->ekin_n = ir->opts.ngtc;
1218     snew(ekinstate->ekinh, ekinstate->ekin_n);
1219     snew(ekinstate->ekinf, ekinstate->ekin_n);
1220     snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1221     ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1222     ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1223     ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1224     ekinstate->dekindl          = 0;
1225     ekinstate->mvcos            = 0;
1226     ekinstate->hasReadEkinState = false;
1227 }
1228
1229 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind)
1230 {
1231     int i;
1232
1233     for (i = 0; i < ekinstate->ekin_n; i++)
1234     {
1235         copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1236         copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1237         copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1238         ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1239         ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1240         ekinstate->vscale_nhc[i]     = ekind->tcstat[i].vscale_nhc;
1241     }
1242
1243     copy_mat(ekind->ekin, ekinstate->ekin_total);
1244     ekinstate->dekindl = ekind->dekindl;
1245     ekinstate->mvcos   = ekind->cosacc.mvcos;
1246 }
1247
1248 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate)
1249 {
1250     int i, n;
1251
1252     if (MASTER(cr))
1253     {
1254         for (i = 0; i < ekinstate->ekin_n; i++)
1255         {
1256             copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1257             copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1258             copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1259             ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1260             ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1261             ekind->tcstat[i].vscale_nhc     = ekinstate->vscale_nhc[i];
1262         }
1263
1264         copy_mat(ekinstate->ekin_total, ekind->ekin);
1265
1266         ekind->dekindl      = ekinstate->dekindl;
1267         ekind->cosacc.mvcos = ekinstate->mvcos;
1268         n                   = ekinstate->ekin_n;
1269     }
1270
1271     if (PAR(cr))
1272     {
1273         gmx_bcast(sizeof(n), &n, cr->mpi_comm_mygroup);
1274         for (i = 0; i < n; i++)
1275         {
1276             gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh[0][0]),
1277                       ekind->tcstat[i].ekinh[0],
1278                       cr->mpi_comm_mygroup);
1279             gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinf[0][0]),
1280                       ekind->tcstat[i].ekinf[0],
1281                       cr->mpi_comm_mygroup);
1282             gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1283                       ekind->tcstat[i].ekinh_old[0],
1284                       cr->mpi_comm_mygroup);
1285
1286             gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1287                       &(ekind->tcstat[i].ekinscalef_nhc),
1288                       cr->mpi_comm_mygroup);
1289             gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1290                       &(ekind->tcstat[i].ekinscaleh_nhc),
1291                       cr->mpi_comm_mygroup);
1292             gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc), &(ekind->tcstat[i].vscale_nhc), cr->mpi_comm_mygroup);
1293         }
1294         gmx_bcast(DIM * DIM * sizeof(ekind->ekin[0][0]), ekind->ekin[0], cr->mpi_comm_mygroup);
1295
1296         gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr->mpi_comm_mygroup);
1297         gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr->mpi_comm_mygroup);
1298     }
1299 }
1300
1301 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom)
1302 {
1303 #if GMX_HAVE_SIMD_UPDATE
1304     constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1305 #else
1306     constexpr int blockSize = 1;
1307 #endif
1308     int numBlocks = (numAtoms + blockSize - 1) / blockSize;
1309
1310     *startAtom = ((numBlocks * threadIndex) / numThreads) * blockSize;
1311     *endAtom   = ((numBlocks * (threadIndex + 1)) / numThreads) * blockSize;
1312     if (threadIndex == numThreads - 1)
1313     {
1314         *endAtom = numAtoms;
1315     }
1316 }
1317
1318 void Update::Impl::update_sd_second_half(const t_inputrec& inputRecord,
1319                                          int64_t           step,
1320                                          real*             dvdlambda,
1321                                          const t_mdatoms*  md,
1322                                          t_state*          state,
1323                                          const t_commrec*  cr,
1324                                          t_nrnb*           nrnb,
1325                                          gmx_wallcycle_t   wcycle,
1326                                          gmx::Constraints* constr,
1327                                          bool              do_log,
1328                                          bool              do_ene)
1329 {
1330     if (!constr)
1331     {
1332         return;
1333     }
1334     if (inputRecord.eI == eiSD1)
1335     {
1336         int homenr = md->homenr;
1337
1338         /* Cast delta_t from double to real to make the integrators faster.
1339          * The only reason for having delta_t double is to get accurate values
1340          * for t=delta_t*step when step is larger than float precision.
1341          * For integration dt the accuracy of real suffices, since with
1342          * integral += dt*integrand the increment is nearly always (much) smaller
1343          * than the integral (and the integrand has real precision).
1344          */
1345         real dt = inputRecord.delta_t;
1346
1347         wallcycle_start(wcycle, ewcUPDATE);
1348
1349         int nth = gmx_omp_nthreads_get(emntUpdate);
1350
1351 #pragma omp parallel for num_threads(nth) schedule(static)
1352         for (int th = 0; th < nth; th++)
1353         {
1354             try
1355             {
1356                 int start_th, end_th;
1357                 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1358
1359                 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>(
1360                         sd_,
1361                         start_th,
1362                         end_th,
1363                         dt,
1364                         inputRecord.opts.nFreeze,
1365                         md->invmass,
1366                         md->ptype,
1367                         md->cFREEZE,
1368                         md->cTC,
1369                         state->x.rvec_array(),
1370                         xp_.rvec_array(),
1371                         state->v.rvec_array(),
1372                         nullptr,
1373                         step,
1374                         inputRecord.ld_seed,
1375                         DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1376             }
1377             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1378         }
1379         inc_nrnb(nrnb, eNR_UPDATE, homenr);
1380         wallcycle_stop(wcycle, ewcUPDATE);
1381
1382         /* Constrain the coordinates upd->xp for half a time step */
1383         bool computeVirial = false;
1384         constr->apply(do_log,
1385                       do_ene,
1386                       step,
1387                       1,
1388                       0.5,
1389                       state->x.arrayRefWithPadding(),
1390                       xp_.arrayRefWithPadding(),
1391                       ArrayRef<RVec>(),
1392                       state->box,
1393                       state->lambda[efptBONDED],
1394                       dvdlambda,
1395                       state->v.arrayRefWithPadding(),
1396                       computeVirial,
1397                       nullptr,
1398                       ConstraintVariable::Positions);
1399     }
1400 }
1401
1402 void Update::Impl::finish_update(const t_inputrec& inputRecord,
1403                                  const t_mdatoms*  md,
1404                                  t_state*          state,
1405                                  gmx_wallcycle_t   wcycle,
1406                                  const bool        haveConstraints)
1407 {
1408     /* NOTE: Currently we always integrate to a temporary buffer and
1409      * then copy the results back here.
1410      */
1411
1412     wallcycle_start_nocount(wcycle, ewcUPDATE);
1413
1414     const int homenr = md->homenr;
1415     auto      xp     = makeConstArrayRef(xp_).subArray(0, homenr);
1416     auto      x      = makeArrayRef(state->x).subArray(0, homenr);
1417
1418     if (md->havePartiallyFrozenAtoms && haveConstraints)
1419     {
1420         /* We have atoms that are frozen along some, but not all dimensions,
1421          * then constraints will have moved them also along the frozen dimensions.
1422          * To freeze such degrees of freedom we do not copy them back here.
1423          */
1424         const ivec* nFreeze = inputRecord.opts.nFreeze;
1425
1426         for (int i = 0; i < homenr; i++)
1427         {
1428             const int g = md->cFREEZE[i];
1429
1430             for (int d = 0; d < DIM; d++)
1431             {
1432                 if (nFreeze[g][d] == 0)
1433                 {
1434                     x[i][d] = xp[i][d];
1435                 }
1436             }
1437         }
1438     }
1439     else
1440     {
1441         /* We have no frozen atoms or fully frozen atoms which have not
1442          * been moved by the update, so we can simply copy all coordinates.
1443          */
1444         int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1445 #pragma omp parallel for num_threads(nth) schedule(static)
1446         for (int i = 0; i < homenr; i++)
1447         {
1448             // Trivial statement, does not throw
1449             x[i] = xp[i];
1450         }
1451     }
1452
1453     wallcycle_stop(wcycle, ewcUPDATE);
1454 }
1455
1456 void Update::Impl::update_coords(const t_inputrec&                                inputRecord,
1457                                  int64_t                                          step,
1458                                  const t_mdatoms*                                 md,
1459                                  t_state*                                         state,
1460                                  const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1461                                  const t_fcdata&                                  fcdata,
1462                                  const gmx_ekindata_t*                            ekind,
1463                                  const matrix                                     M,
1464                                  int                                              updatePart,
1465                                  const t_commrec*                                 cr,
1466                                  const bool                                       haveConstraints)
1467 {
1468     /* Running the velocity half does nothing except for velocity verlet */
1469     if ((updatePart == etrtVELOCITY1 || updatePart == etrtVELOCITY2) && !EI_VV(inputRecord.eI))
1470     {
1471         gmx_incons("update_coords called for velocity without VV integrator");
1472     }
1473
1474     int homenr = md->homenr;
1475
1476     /* Cast to real for faster code, no loss in precision (see comment above) */
1477     real dt = inputRecord.delta_t;
1478
1479     /* We need to update the NMR restraint history when time averaging is used */
1480     if (state->flags & (1 << estDISRE_RM3TAV))
1481     {
1482         update_disres_history(*fcdata.disres, &state->hist);
1483     }
1484     if (state->flags & (1 << estORIRE_DTAV))
1485     {
1486         update_orires_history(*fcdata.orires, &state->hist);
1487     }
1488
1489     /* ############# START The update of velocities and positions ######### */
1490     int nth = gmx_omp_nthreads_get(emntUpdate);
1491
1492 #pragma omp parallel for num_threads(nth) schedule(static)
1493     for (int th = 0; th < nth; th++)
1494     {
1495         try
1496         {
1497             int start_th, end_th;
1498             getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1499
1500             const rvec* x_rvec  = state->x.rvec_array();
1501             rvec*       xp_rvec = xp_.rvec_array();
1502             rvec*       v_rvec  = state->v.rvec_array();
1503             const rvec* f_rvec  = as_rvec_array(f.unpaddedConstArrayRef().data());
1504
1505             switch (inputRecord.eI)
1506             {
1507                 case (eiMD):
1508                     do_update_md(start_th,
1509                                  end_th,
1510                                  dt,
1511                                  step,
1512                                  x_rvec,
1513                                  xp_rvec,
1514                                  v_rvec,
1515                                  f_rvec,
1516                                  inputRecord.etc,
1517                                  inputRecord.epc,
1518                                  inputRecord.nsttcouple,
1519                                  inputRecord.nstpcouple,
1520                                  md,
1521                                  ekind,
1522                                  state->box,
1523                                  state->nosehoover_vxi.data(),
1524                                  M);
1525                     break;
1526                 case (eiSD1):
1527                     do_update_sd(start_th,
1528                                  end_th,
1529                                  dt,
1530                                  step,
1531                                  x_rvec,
1532                                  xp_rvec,
1533                                  v_rvec,
1534                                  f_rvec,
1535                                  inputRecord.opts.nFreeze,
1536                                  md->invmass,
1537                                  md->ptype,
1538                                  md->cFREEZE,
1539                                  md->cTC,
1540                                  inputRecord.ld_seed,
1541                                  cr,
1542                                  sd_,
1543                                  haveConstraints);
1544                     break;
1545                 case (eiBD):
1546                     do_update_bd(start_th,
1547                                  end_th,
1548                                  dt,
1549                                  step,
1550                                  x_rvec,
1551                                  xp_rvec,
1552                                  v_rvec,
1553                                  f_rvec,
1554                                  inputRecord.opts.nFreeze,
1555                                  md->invmass,
1556                                  md->ptype,
1557                                  md->cFREEZE,
1558                                  md->cTC,
1559                                  inputRecord.bd_fric,
1560                                  sd_.bd_rf.data(),
1561                                  inputRecord.ld_seed,
1562                                  DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1563                     break;
1564                 case (eiVV):
1565                 case (eiVVAK):
1566                 {
1567                     gmx_bool bExtended = (inputRecord.etc == etcNOSEHOOVER || inputRecord.epc == epcPARRINELLORAHMAN
1568                                           || inputRecord.epc == epcMTTK);
1569
1570                     /* assuming barostat coupled to group 0 */
1571                     real alpha = 1.0 + DIM / static_cast<real>(inputRecord.opts.nrdf[0]);
1572                     switch (updatePart)
1573                     {
1574                         case etrtVELOCITY1:
1575                         case etrtVELOCITY2:
1576                             do_update_vv_vel(start_th,
1577                                              end_th,
1578                                              dt,
1579                                              inputRecord.opts.nFreeze,
1580                                              md->invmass,
1581                                              md->ptype,
1582                                              md->cFREEZE,
1583                                              v_rvec,
1584                                              f_rvec,
1585                                              bExtended,
1586                                              state->veta,
1587                                              alpha);
1588                             break;
1589                         case etrtPOSITION:
1590                             do_update_vv_pos(start_th,
1591                                              end_th,
1592                                              dt,
1593                                              inputRecord.opts.nFreeze,
1594                                              md->ptype,
1595                                              md->cFREEZE,
1596                                              x_rvec,
1597                                              xp_rvec,
1598                                              v_rvec,
1599                                              bExtended,
1600                                              state->veta);
1601                             break;
1602                     }
1603                     break;
1604                 }
1605                 default: gmx_fatal(FARGS, "Don't know how to update coordinates");
1606             }
1607         }
1608         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1609     }
1610 }
1611
1612 void Update::Impl::update_for_constraint_virial(const t_inputrec& inputRecord,
1613                                                 const t_mdatoms&  md,
1614                                                 const t_state&    state,
1615                                                 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1616                                                 const gmx_ekindata_t& ekind)
1617 {
1618     GMX_ASSERT(inputRecord.eI == eiMD || inputRecord.eI == eiSD1,
1619                "Only leap-frog is supported here");
1620
1621     // Cast to real for faster code, no loss in precision
1622     const real dt = inputRecord.delta_t;
1623
1624     const int nth = gmx_omp_nthreads_get(emntUpdate);
1625
1626 #pragma omp parallel for num_threads(nth) schedule(static)
1627     for (int th = 0; th < nth; th++)
1628     {
1629         try
1630         {
1631             int start_th, end_th;
1632             getThreadAtomRange(nth, th, md.homenr, &start_th, &end_th);
1633
1634             const rvec* x_rvec  = state.x.rvec_array();
1635             rvec*       xp_rvec = xp_.rvec_array();
1636             rvec*       v_rvec  = const_cast<rvec*>(state.v.rvec_array());
1637             const rvec* f_rvec  = as_rvec_array(f.unpaddedConstArrayRef().data());
1638
1639             doUpdateMDDoNotUpdateVelocities(
1640                     start_th, end_th, dt, x_rvec, xp_rvec, v_rvec, f_rvec, md, ekind);
1641         }
1642         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1643     }
1644 }