Make temperature and pressure coupling enums enum classes
[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         switch (accelerationType)
561         {
562             case AccelerationType::none: copy_rvec(v[n], vRel); break;
563             case AccelerationType::cosine:
564                 cosineZ = std::cos(x[n][ZZ] * omega_Z);
565                 vCosine = cosineZ * ekind->cosacc.vcos;
566                 /* Avoid scaling the cosine profile velocity */
567                 copy_rvec(v[n], vRel);
568                 vRel[XX] -= vCosine;
569                 break;
570         }
571
572         real factorNH = 0.0;
573         if (doNoseHoover)
574         {
575             /* Here we account for multiple time stepping, by increasing
576              * the Nose-Hoover correction by nsttcouple
577              * TODO: This can be pre-computed.
578              */
579             factorNH = 0.5 * nsttcouple * dt * nh_vxi[gt];
580         }
581
582         for (int d = 0; d < DIM; d++)
583         {
584             real vNew = (lg * vRel[d]
585                          + (f[n][d] * invMassPerDim[n][d] * dt - factorNH * vRel[d]
586                             - dtPressureCouple * iprod(M[d], vRel)))
587                         / (1 + factorNH);
588             switch (accelerationType)
589             {
590                 case AccelerationType::none: break;
591                 case AccelerationType::cosine:
592                     if (d == XX)
593                     {
594                         /* Add back the mean velocity and apply acceleration */
595                         vNew += vCosine + cosineZ * ekind->cosacc.cos_accel * dt;
596                     }
597                     break;
598             }
599             v[n][d]      = vNew;
600             xprime[n][d] = x[n][d] + vNew * dt;
601         }
602     }
603 }
604
605 /*! \brief Handles the Leap-frog MD x and v integration */
606 static void do_update_md(int         start,
607                          int         nrend,
608                          real        dt,
609                          int64_t     step,
610                          const rvec* gmx_restrict x,
611                          rvec* gmx_restrict xprime,
612                          rvec* gmx_restrict v,
613                          const rvec* gmx_restrict  f,
614                          const TemperatureCoupling etc,
615                          const PressureCoupling    epc,
616                          const int                 nsttcouple,
617                          const int                 nstpcouple,
618                          const t_mdatoms*          md,
619                          const gmx_ekindata_t*     ekind,
620                          const matrix              box,
621                          const double* gmx_restrict nh_vxi,
622                          const matrix               M)
623 {
624     GMX_ASSERT(nrend == start || xprime != x,
625                "For SIMD optimization certain compilers need to have xprime != x");
626
627     /* Note: Berendsen pressure scaling is handled after do_update_md() */
628     bool doTempCouple =
629             (etc != TemperatureCoupling::No && do_per_step(step + nsttcouple - 1, nsttcouple));
630     bool doNoseHoover       = (etc == TemperatureCoupling::NoseHoover && doTempCouple);
631     bool doParrinelloRahman = (epc == PressureCoupling::ParrinelloRahman
632                                && do_per_step(step + nstpcouple - 1, nstpcouple));
633     bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
634
635     real dtPressureCouple = (doParrinelloRahman ? nstpcouple * dt : 0);
636
637     /* NEMD cosine acceleration is applied in updateMDLeapFrogGeneral */
638     bool doAcceleration = (ekind->cosacc.cos_accel != 0);
639
640     if (doNoseHoover || doPROffDiagonal || doAcceleration)
641     {
642         matrix stepM;
643         if (!doParrinelloRahman)
644         {
645             /* We should not apply PR scaling at this step */
646             clear_mat(stepM);
647         }
648         else
649         {
650             copy_mat(M, stepM);
651         }
652
653         if (!doAcceleration)
654         {
655             updateMDLeapfrogGeneral<AccelerationType::none>(
656                     start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
657         }
658         else
659         {
660             updateMDLeapfrogGeneral<AccelerationType::cosine>(
661                     start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
662         }
663     }
664     else
665     {
666         /* Use a simple and thus more efficient integration loop. */
667         /* The simple loop does not check for particle type (so it can
668          * be vectorized), which means we need to clear the velocities
669          * of virtual sites in advance, when present. Note that vsite
670          * velocities are computed after update and constraints from
671          * their displacement.
672          */
673         if (md->haveVsites)
674         {
675             /* Note: The overhead of this loop is completely neligible */
676             clearVsiteVelocities(start, nrend, md->ptype, v);
677         }
678
679         /* We determine if we have a single T-coupling lambda value for all
680          * atoms. That allows for better SIMD acceleration in the template.
681          * If we do not do temperature coupling (in the run or this step),
682          * all scaling values are 1, so we effectively have a single value.
683          */
684         bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
685
686         /* Extract some pointers needed by all cases */
687         const unsigned short*             cTC           = md->cTC;
688         gmx::ArrayRef<const t_grp_tcstat> tcstat        = ekind->tcstat;
689         const rvec*                       invMassPerDim = md->invMassPerDim;
690
691         if (doParrinelloRahman)
692         {
693             GMX_ASSERT(!doPROffDiagonal,
694                        "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling "
695                        "matrices");
696
697             rvec diagM;
698             for (int d = 0; d < DIM; d++)
699             {
700                 diagM[d] = M[d][d];
701             }
702
703             if (haveSingleTempScaleValue)
704             {
705                 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::diagonal>(
706                         start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
707             }
708             else
709             {
710                 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::diagonal>(
711                         start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
712             }
713         }
714         else
715         {
716             if (haveSingleTempScaleValue)
717             {
718                 /* Note that modern compilers are pretty good at vectorizing
719                  * updateMDLeapfrogSimple(). But the SIMD version will still
720                  * be faster because invMass lowers the cache pressure
721                  * compared to invMassPerDim.
722                  */
723 #if GMX_HAVE_SIMD_UPDATE
724                 /* Check if we can use invmass instead of invMassPerDim */
725                 if (!md->havePartiallyFrozenAtoms)
726                 {
727                     updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::yes>(
728                             start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
729                 }
730                 else
731 #endif
732                 {
733                     updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
734                             start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
735                 }
736             }
737             else
738             {
739                 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::no>(
740                         start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
741             }
742         }
743     }
744 }
745 /*! \brief Handles the Leap-frog MD x and v integration */
746 static void doUpdateMDDoNotUpdateVelocities(int         start,
747                                             int         nrend,
748                                             real        dt,
749                                             const rvec* gmx_restrict x,
750                                             rvec* gmx_restrict xprime,
751                                             rvec* gmx_restrict v,
752                                             const rvec* gmx_restrict f,
753                                             const t_mdatoms&         md,
754                                             const gmx_ekindata_t&    ekind)
755 {
756     GMX_ASSERT(nrend == start || xprime != x,
757                "For SIMD optimization certain compilers need to have xprime != x");
758
759     gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind.tcstat;
760
761     /* Check if we can use invmass instead of invMassPerDim */
762 #if GMX_HAVE_SIMD_UPDATE
763     if (!md.havePartiallyFrozenAtoms)
764     {
765         updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::no>(
766                 start, nrend, dt, md.invmass, tcstat, x, xprime, v, f);
767     }
768     else
769 #endif
770     {
771         updateMDLeapfrogSimple<StoreUpdatedVelocities::no, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
772                 start, nrend, dt, dt, md.invMassPerDim, tcstat, nullptr, nullptr, x, xprime, v, f);
773     }
774 }
775
776 static void do_update_vv_vel(int                  start,
777                              int                  nrend,
778                              real                 dt,
779                              const ivec           nFreeze[],
780                              const real           invmass[],
781                              const unsigned short ptype[],
782                              const unsigned short cFREEZE[],
783                              rvec                 v[],
784                              const rvec           f[],
785                              gmx_bool             bExtended,
786                              real                 veta,
787                              real                 alpha)
788 {
789     int  gf = 0;
790     int  n, d;
791     real g, mv1, mv2;
792
793     if (bExtended)
794     {
795         g   = 0.25 * dt * veta * alpha;
796         mv1 = std::exp(-g);
797         mv2 = gmx::series_sinhx(g);
798     }
799     else
800     {
801         mv1 = 1.0;
802         mv2 = 1.0;
803     }
804     for (n = start; n < nrend; n++)
805     {
806         real w_dt = invmass[n] * dt;
807         if (cFREEZE)
808         {
809             gf = cFREEZE[n];
810         }
811
812         for (d = 0; d < DIM; d++)
813         {
814             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
815             {
816                 v[n][d] = mv1 * (mv1 * v[n][d] + 0.5 * (w_dt * mv2 * f[n][d]));
817             }
818             else
819             {
820                 v[n][d] = 0.0;
821             }
822         }
823     }
824 } /* do_update_vv_vel */
825
826 static void do_update_vv_pos(int                  start,
827                              int                  nrend,
828                              real                 dt,
829                              const ivec           nFreeze[],
830                              const unsigned short ptype[],
831                              const unsigned short cFREEZE[],
832                              const rvec           x[],
833                              rvec                 xprime[],
834                              const rvec           v[],
835                              gmx_bool             bExtended,
836                              real                 veta)
837 {
838     int  gf = 0;
839     int  n, d;
840     real g, mr1, mr2;
841
842     /* Would it make more sense if Parrinello-Rahman was put here? */
843     if (bExtended)
844     {
845         g   = 0.5 * dt * veta;
846         mr1 = std::exp(g);
847         mr2 = gmx::series_sinhx(g);
848     }
849     else
850     {
851         mr1 = 1.0;
852         mr2 = 1.0;
853     }
854
855     for (n = start; n < nrend; n++)
856     {
857
858         if (cFREEZE)
859         {
860             gf = cFREEZE[n];
861         }
862
863         for (d = 0; d < DIM; d++)
864         {
865             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
866             {
867                 xprime[n][d] = mr1 * (mr1 * x[n][d] + mr2 * dt * v[n][d]);
868             }
869             else
870             {
871                 xprime[n][d] = x[n][d];
872             }
873         }
874     }
875 } /* do_update_vv_pos */
876
877 gmx_stochd_t::gmx_stochd_t(const t_inputrec& inputRecord)
878 {
879     const t_grpopts* opts = &inputRecord.opts;
880     const int        ngtc = opts->ngtc;
881
882     if (inputRecord.eI == eiBD)
883     {
884         bd_rf.resize(ngtc);
885     }
886     else if (EI_SD(inputRecord.eI))
887     {
888         sdc.resize(ngtc);
889         sdsig.resize(ngtc);
890
891         for (int gt = 0; gt < ngtc; gt++)
892         {
893             if (opts->tau_t[gt] > 0)
894             {
895                 sdc[gt].em = std::exp(-inputRecord.delta_t / opts->tau_t[gt]);
896             }
897             else
898             {
899                 /* No friction and noise on this group */
900                 sdc[gt].em = 1;
901             }
902         }
903     }
904     else if (ETC_ANDERSEN(inputRecord.etc))
905     {
906         randomize_group.resize(ngtc);
907         boltzfac.resize(ngtc);
908
909         /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
910         /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
911
912         for (int gt = 0; gt < ngtc; gt++)
913         {
914             real reft = std::max<real>(0, opts->ref_t[gt]);
915             if ((opts->tau_t[gt] > 0)
916                 && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
917             {
918                 randomize_group[gt] = true;
919                 boltzfac[gt]        = BOLTZ * opts->ref_t[gt];
920             }
921             else
922             {
923                 randomize_group[gt] = false;
924             }
925         }
926     }
927 }
928
929 void Update::Impl::update_temperature_constants(const t_inputrec& inputRecord)
930 {
931     if (inputRecord.eI == eiBD)
932     {
933         if (inputRecord.bd_fric != 0)
934         {
935             for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
936             {
937                 sd_.bd_rf[gt] = std::sqrt(2.0 * BOLTZ * inputRecord.opts.ref_t[gt]
938                                           / (inputRecord.bd_fric * inputRecord.delta_t));
939             }
940         }
941         else
942         {
943             for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
944             {
945                 sd_.bd_rf[gt] = std::sqrt(2.0 * BOLTZ * inputRecord.opts.ref_t[gt]);
946             }
947         }
948     }
949     if (inputRecord.eI == eiSD1)
950     {
951         for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
952         {
953             real kT = BOLTZ * inputRecord.opts.ref_t[gt];
954             /* The mass is accounted for later, since this differs per atom */
955             sd_.sdsig[gt].V = std::sqrt(kT * (1 - sd_.sdc[gt].em * sd_.sdc[gt].em));
956         }
957     }
958 }
959
960 Update::Impl::Impl(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
961     sd_(inputRecord),
962     deform_(boxDeformation)
963 {
964     update_temperature_constants(inputRecord);
965     xp_.resizeWithPadding(0);
966 }
967
968 void Update::setNumAtoms(int numAtoms)
969 {
970
971     impl_->xp()->resizeWithPadding(numAtoms);
972 }
973
974 /*! \brief Sets the SD update type */
975 enum class SDUpdate : int
976 {
977     ForcesOnly,
978     FrictionAndNoiseOnly,
979     Combined
980 };
981
982 /*! \brief SD integrator update
983  *
984  * Two phases are required in the general case of a constrained
985  * update, the first phase from the contribution of forces, before
986  * applying constraints, and then a second phase applying the friction
987  * and noise, and then further constraining. For details, see
988  * Goga2012.
989  *
990  * Without constraints, the two phases can be combined, for
991  * efficiency.
992  *
993  * Thus three instantiations of this templated function will be made,
994  * two with only one contribution, and one with both contributions. */
995 template<SDUpdate updateType>
996 static void doSDUpdateGeneral(const gmx_stochd_t&  sd,
997                               int                  start,
998                               int                  nrend,
999                               real                 dt,
1000                               const ivec           nFreeze[],
1001                               const real           invmass[],
1002                               const unsigned short ptype[],
1003                               const unsigned short cFREEZE[],
1004                               const unsigned short cTC[],
1005                               const rvec           x[],
1006                               rvec                 xprime[],
1007                               rvec                 v[],
1008                               const rvec           f[],
1009                               int64_t              step,
1010                               int                  seed,
1011                               const int*           gatindex)
1012 {
1013     // cTC and cFREEZE can be nullptr any time, but various
1014     // instantiations do not make sense with particular pointer
1015     // values.
1016     if (updateType == SDUpdate::ForcesOnly)
1017     {
1018         GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
1019         GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
1020     }
1021     if (updateType == SDUpdate::FrictionAndNoiseOnly)
1022     {
1023         GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
1024     }
1025     if (updateType == SDUpdate::Combined)
1026     {
1027         GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
1028     }
1029
1030     // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
1031     gmx::ThreeFry2x64<0>                       rng(seed, gmx::RandomDomain::UpdateCoordinates);
1032     gmx::TabulatedNormalDistribution<real, 14> dist;
1033
1034     for (int n = start; n < nrend; n++)
1035     {
1036         int globalAtomIndex = gatindex ? gatindex[n] : n;
1037         rng.restart(step, globalAtomIndex);
1038         dist.reset();
1039
1040         real inverseMass = invmass[n];
1041         real invsqrtMass = std::sqrt(inverseMass);
1042
1043         int freezeGroup      = cFREEZE ? cFREEZE[n] : 0;
1044         int temperatureGroup = cTC ? cTC[n] : 0;
1045
1046         for (int d = 0; d < DIM; d++)
1047         {
1048             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[freezeGroup][d])
1049             {
1050                 if (updateType == SDUpdate::ForcesOnly)
1051                 {
1052                     real vn = v[n][d] + inverseMass * f[n][d] * dt;
1053                     v[n][d] = vn;
1054                     // Simple position update.
1055                     xprime[n][d] = x[n][d] + v[n][d] * dt;
1056                 }
1057                 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
1058                 {
1059                     real vn = v[n][d];
1060                     v[n][d] = (vn * sd.sdc[temperatureGroup].em
1061                                + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1062                     // The previous phase already updated the
1063                     // positions with a full v*dt term that must
1064                     // now be half removed.
1065                     xprime[n][d] = xprime[n][d] + 0.5 * (v[n][d] - vn) * dt;
1066                 }
1067                 else
1068                 {
1069                     real vn = v[n][d] + inverseMass * f[n][d] * dt;
1070                     v[n][d] = (vn * sd.sdc[temperatureGroup].em
1071                                + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1072                     // Here we include half of the friction+noise
1073                     // update of v into the position update.
1074                     xprime[n][d] = x[n][d] + 0.5 * (vn + v[n][d]) * dt;
1075                 }
1076             }
1077             else
1078             {
1079                 // When using constraints, the update is split into
1080                 // two phases, but we only need to zero the update of
1081                 // virtual, shell or frozen particles in at most one
1082                 // of the phases.
1083                 if (updateType != SDUpdate::FrictionAndNoiseOnly)
1084                 {
1085                     v[n][d]      = 0.0;
1086                     xprime[n][d] = x[n][d];
1087                 }
1088             }
1089         }
1090     }
1091 }
1092
1093 static void do_update_sd(int         start,
1094                          int         nrend,
1095                          real        dt,
1096                          int64_t     step,
1097                          const rvec* gmx_restrict x,
1098                          rvec* gmx_restrict xprime,
1099                          rvec* gmx_restrict v,
1100                          const rvec* gmx_restrict f,
1101                          const ivec               nFreeze[],
1102                          const real               invmass[],
1103                          const unsigned short     ptype[],
1104                          const unsigned short     cFREEZE[],
1105                          const unsigned short     cTC[],
1106                          int                      seed,
1107                          const t_commrec*         cr,
1108                          const gmx_stochd_t&      sd,
1109                          bool                     haveConstraints)
1110 {
1111     if (haveConstraints)
1112     {
1113         // With constraints, the SD update is done in 2 parts
1114         doSDUpdateGeneral<SDUpdate::ForcesOnly>(
1115                 sd, start, nrend, dt, nFreeze, invmass, ptype, cFREEZE, nullptr, x, xprime, v, f, step, seed, nullptr);
1116     }
1117     else
1118     {
1119         doSDUpdateGeneral<SDUpdate::Combined>(sd,
1120                                               start,
1121                                               nrend,
1122                                               dt,
1123                                               nFreeze,
1124                                               invmass,
1125                                               ptype,
1126                                               cFREEZE,
1127                                               cTC,
1128                                               x,
1129                                               xprime,
1130                                               v,
1131                                               f,
1132                                               step,
1133                                               seed,
1134                                               DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1135     }
1136 }
1137
1138 static void do_update_bd(int         start,
1139                          int         nrend,
1140                          real        dt,
1141                          int64_t     step,
1142                          const rvec* gmx_restrict x,
1143                          rvec* gmx_restrict xprime,
1144                          rvec* gmx_restrict v,
1145                          const rvec* gmx_restrict f,
1146                          const ivec               nFreeze[],
1147                          const real               invmass[],
1148                          const unsigned short     ptype[],
1149                          const unsigned short     cFREEZE[],
1150                          const unsigned short     cTC[],
1151                          real                     friction_coefficient,
1152                          const real*              rf,
1153                          int                      seed,
1154                          const int*               gatindex)
1155 {
1156     /* note -- these appear to be full step velocities . . .  */
1157     int  gf = 0, gt = 0;
1158     real vn;
1159     real invfr = 0;
1160     int  n, d;
1161     // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1162     // Each 64-bit value is enough for 4 normal distribution table numbers.
1163     gmx::ThreeFry2x64<0>                       rng(seed, gmx::RandomDomain::UpdateCoordinates);
1164     gmx::TabulatedNormalDistribution<real, 14> dist;
1165
1166     if (friction_coefficient != 0)
1167     {
1168         invfr = 1.0 / friction_coefficient;
1169     }
1170
1171     for (n = start; (n < nrend); n++)
1172     {
1173         int ng = gatindex ? gatindex[n] : n;
1174
1175         rng.restart(step, ng);
1176         dist.reset();
1177
1178         if (cFREEZE)
1179         {
1180             gf = cFREEZE[n];
1181         }
1182         if (cTC)
1183         {
1184             gt = cTC[n];
1185         }
1186         for (d = 0; (d < DIM); d++)
1187         {
1188             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
1189             {
1190                 if (friction_coefficient != 0)
1191                 {
1192                     vn = invfr * f[n][d] + rf[gt] * dist(rng);
1193                 }
1194                 else
1195                 {
1196                     /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1197                     vn = 0.5 * invmass[n] * f[n][d] * dt
1198                          + std::sqrt(0.5 * invmass[n]) * rf[gt] * dist(rng);
1199                 }
1200
1201                 v[n][d]      = vn;
1202                 xprime[n][d] = x[n][d] + vn * dt;
1203             }
1204             else
1205             {
1206                 v[n][d]      = 0.0;
1207                 xprime[n][d] = x[n][d];
1208             }
1209         }
1210     }
1211 }
1212
1213 extern void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir)
1214 {
1215     ekinstate->ekin_n = ir->opts.ngtc;
1216     snew(ekinstate->ekinh, ekinstate->ekin_n);
1217     snew(ekinstate->ekinf, ekinstate->ekin_n);
1218     snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1219     ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1220     ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1221     ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1222     ekinstate->dekindl          = 0;
1223     ekinstate->mvcos            = 0;
1224     ekinstate->hasReadEkinState = false;
1225 }
1226
1227 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind)
1228 {
1229     int i;
1230
1231     for (i = 0; i < ekinstate->ekin_n; i++)
1232     {
1233         copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1234         copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1235         copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1236         ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1237         ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1238         ekinstate->vscale_nhc[i]     = ekind->tcstat[i].vscale_nhc;
1239     }
1240
1241     copy_mat(ekind->ekin, ekinstate->ekin_total);
1242     ekinstate->dekindl = ekind->dekindl;
1243     ekinstate->mvcos   = ekind->cosacc.mvcos;
1244 }
1245
1246 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate)
1247 {
1248     int i, n;
1249
1250     if (MASTER(cr))
1251     {
1252         for (i = 0; i < ekinstate->ekin_n; i++)
1253         {
1254             copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1255             copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1256             copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1257             ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1258             ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1259             ekind->tcstat[i].vscale_nhc     = ekinstate->vscale_nhc[i];
1260         }
1261
1262         copy_mat(ekinstate->ekin_total, ekind->ekin);
1263
1264         ekind->dekindl      = ekinstate->dekindl;
1265         ekind->cosacc.mvcos = ekinstate->mvcos;
1266         n                   = ekinstate->ekin_n;
1267     }
1268
1269     if (PAR(cr))
1270     {
1271         gmx_bcast(sizeof(n), &n, cr->mpi_comm_mygroup);
1272         for (i = 0; i < n; i++)
1273         {
1274             gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh[0][0]),
1275                       ekind->tcstat[i].ekinh[0],
1276                       cr->mpi_comm_mygroup);
1277             gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinf[0][0]),
1278                       ekind->tcstat[i].ekinf[0],
1279                       cr->mpi_comm_mygroup);
1280             gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1281                       ekind->tcstat[i].ekinh_old[0],
1282                       cr->mpi_comm_mygroup);
1283
1284             gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1285                       &(ekind->tcstat[i].ekinscalef_nhc),
1286                       cr->mpi_comm_mygroup);
1287             gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1288                       &(ekind->tcstat[i].ekinscaleh_nhc),
1289                       cr->mpi_comm_mygroup);
1290             gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc), &(ekind->tcstat[i].vscale_nhc), cr->mpi_comm_mygroup);
1291         }
1292         gmx_bcast(DIM * DIM * sizeof(ekind->ekin[0][0]), ekind->ekin[0], cr->mpi_comm_mygroup);
1293
1294         gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr->mpi_comm_mygroup);
1295         gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr->mpi_comm_mygroup);
1296     }
1297 }
1298
1299 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom)
1300 {
1301 #if GMX_HAVE_SIMD_UPDATE
1302     constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1303 #else
1304     constexpr int blockSize = 1;
1305 #endif
1306     int numBlocks = (numAtoms + blockSize - 1) / blockSize;
1307
1308     *startAtom = ((numBlocks * threadIndex) / numThreads) * blockSize;
1309     *endAtom   = ((numBlocks * (threadIndex + 1)) / numThreads) * blockSize;
1310     if (threadIndex == numThreads - 1)
1311     {
1312         *endAtom = numAtoms;
1313     }
1314 }
1315
1316 void Update::Impl::update_sd_second_half(const t_inputrec& inputRecord,
1317                                          int64_t           step,
1318                                          real*             dvdlambda,
1319                                          const t_mdatoms*  md,
1320                                          t_state*          state,
1321                                          const t_commrec*  cr,
1322                                          t_nrnb*           nrnb,
1323                                          gmx_wallcycle_t   wcycle,
1324                                          gmx::Constraints* constr,
1325                                          bool              do_log,
1326                                          bool              do_ene)
1327 {
1328     if (!constr)
1329     {
1330         return;
1331     }
1332     if (inputRecord.eI == eiSD1)
1333     {
1334         int homenr = md->homenr;
1335
1336         /* Cast delta_t from double to real to make the integrators faster.
1337          * The only reason for having delta_t double is to get accurate values
1338          * for t=delta_t*step when step is larger than float precision.
1339          * For integration dt the accuracy of real suffices, since with
1340          * integral += dt*integrand the increment is nearly always (much) smaller
1341          * than the integral (and the integrand has real precision).
1342          */
1343         real dt = inputRecord.delta_t;
1344
1345         wallcycle_start(wcycle, ewcUPDATE);
1346
1347         int nth = gmx_omp_nthreads_get(emntUpdate);
1348
1349 #pragma omp parallel for num_threads(nth) schedule(static)
1350         for (int th = 0; th < nth; th++)
1351         {
1352             try
1353             {
1354                 int start_th, end_th;
1355                 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1356
1357                 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>(
1358                         sd_,
1359                         start_th,
1360                         end_th,
1361                         dt,
1362                         inputRecord.opts.nFreeze,
1363                         md->invmass,
1364                         md->ptype,
1365                         md->cFREEZE,
1366                         md->cTC,
1367                         state->x.rvec_array(),
1368                         xp_.rvec_array(),
1369                         state->v.rvec_array(),
1370                         nullptr,
1371                         step,
1372                         inputRecord.ld_seed,
1373                         DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1374             }
1375             GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1376         }
1377         inc_nrnb(nrnb, eNR_UPDATE, homenr);
1378         wallcycle_stop(wcycle, ewcUPDATE);
1379
1380         /* Constrain the coordinates upd->xp for half a time step */
1381         bool computeVirial = false;
1382         constr->apply(do_log,
1383                       do_ene,
1384                       step,
1385                       1,
1386                       0.5,
1387                       state->x.arrayRefWithPadding(),
1388                       xp_.arrayRefWithPadding(),
1389                       ArrayRef<RVec>(),
1390                       state->box,
1391                       state->lambda[efptBONDED],
1392                       dvdlambda,
1393                       state->v.arrayRefWithPadding(),
1394                       computeVirial,
1395                       nullptr,
1396                       ConstraintVariable::Positions);
1397     }
1398 }
1399
1400 void Update::Impl::finish_update(const t_inputrec& inputRecord,
1401                                  const t_mdatoms*  md,
1402                                  t_state*          state,
1403                                  gmx_wallcycle_t   wcycle,
1404                                  const bool        haveConstraints)
1405 {
1406     /* NOTE: Currently we always integrate to a temporary buffer and
1407      * then copy the results back here.
1408      */
1409
1410     wallcycle_start_nocount(wcycle, ewcUPDATE);
1411
1412     const int homenr = md->homenr;
1413     auto      xp     = makeConstArrayRef(xp_).subArray(0, homenr);
1414     auto      x      = makeArrayRef(state->x).subArray(0, homenr);
1415
1416     if (md->havePartiallyFrozenAtoms && haveConstraints)
1417     {
1418         /* We have atoms that are frozen along some, but not all dimensions,
1419          * then constraints will have moved them also along the frozen dimensions.
1420          * To freeze such degrees of freedom we do not copy them back here.
1421          */
1422         const ivec* nFreeze = inputRecord.opts.nFreeze;
1423
1424         for (int i = 0; i < homenr; i++)
1425         {
1426             const int g = md->cFREEZE[i];
1427
1428             for (int d = 0; d < DIM; d++)
1429             {
1430                 if (nFreeze[g][d] == 0)
1431                 {
1432                     x[i][d] = xp[i][d];
1433                 }
1434             }
1435         }
1436     }
1437     else
1438     {
1439         /* We have no frozen atoms or fully frozen atoms which have not
1440          * been moved by the update, so we can simply copy all coordinates.
1441          */
1442         int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1443 #pragma omp parallel for num_threads(nth) schedule(static)
1444         for (int i = 0; i < homenr; i++)
1445         {
1446             // Trivial statement, does not throw
1447             x[i] = xp[i];
1448         }
1449     }
1450
1451     wallcycle_stop(wcycle, ewcUPDATE);
1452 }
1453
1454 void Update::Impl::update_coords(const t_inputrec&                                inputRecord,
1455                                  int64_t                                          step,
1456                                  const t_mdatoms*                                 md,
1457                                  t_state*                                         state,
1458                                  const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1459                                  const t_fcdata&                                  fcdata,
1460                                  const gmx_ekindata_t*                            ekind,
1461                                  const matrix                                     M,
1462                                  int                                              updatePart,
1463                                  const t_commrec*                                 cr,
1464                                  const bool                                       haveConstraints)
1465 {
1466     /* Running the velocity half does nothing except for velocity verlet */
1467     if ((updatePart == etrtVELOCITY1 || updatePart == etrtVELOCITY2) && !EI_VV(inputRecord.eI))
1468     {
1469         gmx_incons("update_coords called for velocity without VV integrator");
1470     }
1471
1472     int homenr = md->homenr;
1473
1474     /* Cast to real for faster code, no loss in precision (see comment above) */
1475     real dt = inputRecord.delta_t;
1476
1477     /* We need to update the NMR restraint history when time averaging is used */
1478     if (state->flags & (1 << estDISRE_RM3TAV))
1479     {
1480         update_disres_history(*fcdata.disres, &state->hist);
1481     }
1482     if (state->flags & (1 << estORIRE_DTAV))
1483     {
1484         update_orires_history(*fcdata.orires, &state->hist);
1485     }
1486
1487     /* ############# START The update of velocities and positions ######### */
1488     int nth = gmx_omp_nthreads_get(emntUpdate);
1489
1490 #pragma omp parallel for num_threads(nth) schedule(static)
1491     for (int th = 0; th < nth; th++)
1492     {
1493         try
1494         {
1495             int start_th, end_th;
1496             getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1497
1498             const rvec* x_rvec  = state->x.rvec_array();
1499             rvec*       xp_rvec = xp_.rvec_array();
1500             rvec*       v_rvec  = state->v.rvec_array();
1501             const rvec* f_rvec  = as_rvec_array(f.unpaddedConstArrayRef().data());
1502
1503             switch (inputRecord.eI)
1504             {
1505                 case (eiMD):
1506                     do_update_md(start_th,
1507                                  end_th,
1508                                  dt,
1509                                  step,
1510                                  x_rvec,
1511                                  xp_rvec,
1512                                  v_rvec,
1513                                  f_rvec,
1514                                  inputRecord.etc,
1515                                  inputRecord.epc,
1516                                  inputRecord.nsttcouple,
1517                                  inputRecord.nstpcouple,
1518                                  md,
1519                                  ekind,
1520                                  state->box,
1521                                  state->nosehoover_vxi.data(),
1522                                  M);
1523                     break;
1524                 case (eiSD1):
1525                     do_update_sd(start_th,
1526                                  end_th,
1527                                  dt,
1528                                  step,
1529                                  x_rvec,
1530                                  xp_rvec,
1531                                  v_rvec,
1532                                  f_rvec,
1533                                  inputRecord.opts.nFreeze,
1534                                  md->invmass,
1535                                  md->ptype,
1536                                  md->cFREEZE,
1537                                  md->cTC,
1538                                  inputRecord.ld_seed,
1539                                  cr,
1540                                  sd_,
1541                                  haveConstraints);
1542                     break;
1543                 case (eiBD):
1544                     do_update_bd(start_th,
1545                                  end_th,
1546                                  dt,
1547                                  step,
1548                                  x_rvec,
1549                                  xp_rvec,
1550                                  v_rvec,
1551                                  f_rvec,
1552                                  inputRecord.opts.nFreeze,
1553                                  md->invmass,
1554                                  md->ptype,
1555                                  md->cFREEZE,
1556                                  md->cTC,
1557                                  inputRecord.bd_fric,
1558                                  sd_.bd_rf.data(),
1559                                  inputRecord.ld_seed,
1560                                  DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1561                     break;
1562                 case (eiVV):
1563                 case (eiVVAK):
1564                 {
1565                     gmx_bool bExtended = (inputRecord.etc == TemperatureCoupling::NoseHoover
1566                                           || inputRecord.epc == PressureCoupling::ParrinelloRahman
1567                                           || inputRecord.epc == PressureCoupling::Mttk);
1568
1569                     /* assuming barostat coupled to group 0 */
1570                     real alpha = 1.0 + DIM / static_cast<real>(inputRecord.opts.nrdf[0]);
1571                     switch (updatePart)
1572                     {
1573                         case etrtVELOCITY1:
1574                         case etrtVELOCITY2:
1575                             do_update_vv_vel(start_th,
1576                                              end_th,
1577                                              dt,
1578                                              inputRecord.opts.nFreeze,
1579                                              md->invmass,
1580                                              md->ptype,
1581                                              md->cFREEZE,
1582                                              v_rvec,
1583                                              f_rvec,
1584                                              bExtended,
1585                                              state->veta,
1586                                              alpha);
1587                             break;
1588                         case etrtPOSITION:
1589                             do_update_vv_pos(start_th,
1590                                              end_th,
1591                                              dt,
1592                                              inputRecord.opts.nFreeze,
1593                                              md->ptype,
1594                                              md->cFREEZE,
1595                                              x_rvec,
1596                                              xp_rvec,
1597                                              v_rvec,
1598                                              bExtended,
1599                                              state->veta);
1600                             break;
1601                     }
1602                     break;
1603                 }
1604                 default: gmx_fatal(FARGS, "Don't know how to update coordinates");
1605             }
1606         }
1607         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1608     }
1609 }
1610
1611 void Update::Impl::update_for_constraint_virial(const t_inputrec& inputRecord,
1612                                                 const t_mdatoms&  md,
1613                                                 const t_state&    state,
1614                                                 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1615                                                 const gmx_ekindata_t& ekind)
1616 {
1617     GMX_ASSERT(inputRecord.eI == eiMD || inputRecord.eI == eiSD1,
1618                "Only leap-frog is supported here");
1619
1620     // Cast to real for faster code, no loss in precision
1621     const real dt = inputRecord.delta_t;
1622
1623     const int nth = gmx_omp_nthreads_get(emntUpdate);
1624
1625 #pragma omp parallel for num_threads(nth) schedule(static)
1626     for (int th = 0; th < nth; th++)
1627     {
1628         try
1629         {
1630             int start_th, end_th;
1631             getThreadAtomRange(nth, th, md.homenr, &start_th, &end_th);
1632
1633             const rvec* x_rvec  = state.x.rvec_array();
1634             rvec*       xp_rvec = xp_.rvec_array();
1635             rvec*       v_rvec  = const_cast<rvec*>(state.v.rvec_array());
1636             const rvec* f_rvec  = as_rvec_array(f.unpaddedConstArrayRef().data());
1637
1638             doUpdateMDDoNotUpdateVelocities(
1639                     start_th, end_th, dt, x_rvec, xp_rvec, v_rvec, f_rvec, md, ekind);
1640         }
1641         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1642     }
1643 }