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