Merge release-5-0 into master
[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, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/update.h"
40
41 #include <math.h>
42 #include <stdio.h>
43
44 #include <algorithm>
45
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/legacyheaders/constr.h"
48 #include "gromacs/legacyheaders/disre.h"
49 #include "gromacs/legacyheaders/force.h"
50 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
51 #include "gromacs/legacyheaders/macros.h"
52 #include "gromacs/legacyheaders/mdrun.h"
53 #include "gromacs/legacyheaders/names.h"
54 #include "gromacs/legacyheaders/nrnb.h"
55 #include "gromacs/legacyheaders/orires.h"
56 #include "gromacs/legacyheaders/tgroup.h"
57 #include "gromacs/legacyheaders/txtdump.h"
58 #include "gromacs/legacyheaders/typedefs.h"
59 #include "gromacs/legacyheaders/types/commrec.h"
60 #include "gromacs/math/units.h"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/pbcutil/mshift.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/pulling/pull.h"
65 #include "gromacs/random/random.h"
66 #include "gromacs/timing/wallcycle.h"
67 #include "gromacs/utility/futil.h"
68 #include "gromacs/utility/gmxomp.h"
69 #include "gromacs/utility/smalloc.h"
70
71 /*For debugging, start at v(-dt/2) for velolcity verlet -- uncomment next line */
72 /*#define STARTFROMDT2*/
73
74 typedef struct {
75     double gdt;
76     double eph;
77     double emh;
78     double em;
79     double b;
80     double c;
81     double d;
82 } gmx_sd_const_t;
83
84 typedef struct {
85     real V;
86     real X;
87     real Yv;
88     real Yx;
89 } gmx_sd_sigma_t;
90
91 typedef struct {
92     /* BD stuff */
93     real           *bd_rf;
94     /* SD stuff */
95     gmx_sd_const_t *sdc;
96     gmx_sd_sigma_t *sdsig;
97     rvec           *sd_V;
98     int             sd_V_nalloc;
99     /* andersen temperature control stuff */
100     gmx_bool       *randomize_group;
101     real           *boltzfac;
102 } gmx_stochd_t;
103
104 typedef struct gmx_update
105 {
106     gmx_stochd_t *sd;
107     /* xprime for constraint algorithms */
108     rvec         *xp;
109     int           xp_nalloc;
110
111     /* Variables for the deform algorithm */
112     gmx_int64_t     deformref_step;
113     matrix          deformref_box;
114 } t_gmx_update;
115
116
117 static void do_update_md(int start, int nrend, double dt,
118                          t_grp_tcstat *tcstat,
119                          double nh_vxi[],
120                          gmx_bool bNEMD, t_grp_acc *gstat, rvec accel[],
121                          ivec nFreeze[],
122                          real invmass[],
123                          unsigned short ptype[], unsigned short cFREEZE[],
124                          unsigned short cACC[], unsigned short cTC[],
125                          rvec x[], rvec xprime[], rvec v[],
126                          rvec f[], matrix M,
127                          gmx_bool bNH, gmx_bool bPR)
128 {
129     double imass, w_dt;
130     int    gf = 0, ga = 0, gt = 0;
131     rvec   vrel;
132     real   vn, vv, va, vb, vnrel;
133     real   lg, vxi = 0, u;
134     int    n, d;
135
136     if (bNH || bPR)
137     {
138         /* Update with coupling to extended ensembles, used for
139          * Nose-Hoover and Parrinello-Rahman coupling
140          * Nose-Hoover uses the reversible leap-frog integrator from
141          * Holian et al. Phys Rev E 52(3) : 2338, 1995
142          */
143         for (n = start; n < nrend; n++)
144         {
145             imass = invmass[n];
146             if (cFREEZE)
147             {
148                 gf   = cFREEZE[n];
149             }
150             if (cACC)
151             {
152                 ga   = cACC[n];
153             }
154             if (cTC)
155             {
156                 gt   = cTC[n];
157             }
158             lg   = tcstat[gt].lambda;
159             if (bNH)
160             {
161                 vxi   = nh_vxi[gt];
162             }
163             rvec_sub(v[n], gstat[ga].u, vrel);
164
165             for (d = 0; d < DIM; d++)
166             {
167                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
168                 {
169                     vnrel = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
170                                               - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
171                     /* do not scale the mean velocities u */
172                     vn             = gstat[ga].u[d] + accel[ga][d]*dt + vnrel;
173                     v[n][d]        = vn;
174                     xprime[n][d]   = x[n][d]+vn*dt;
175                 }
176                 else
177                 {
178                     v[n][d]        = 0.0;
179                     xprime[n][d]   = x[n][d];
180                 }
181             }
182         }
183     }
184     else if (cFREEZE != NULL ||
185              nFreeze[0][XX] || nFreeze[0][YY] || nFreeze[0][ZZ] ||
186              bNEMD)
187     {
188         /* Update with Berendsen/v-rescale coupling and freeze or NEMD */
189         for (n = start; n < nrend; n++)
190         {
191             w_dt = invmass[n]*dt;
192             if (cFREEZE)
193             {
194                 gf   = cFREEZE[n];
195             }
196             if (cACC)
197             {
198                 ga   = cACC[n];
199             }
200             if (cTC)
201             {
202                 gt   = cTC[n];
203             }
204             lg   = tcstat[gt].lambda;
205
206             for (d = 0; d < DIM; d++)
207             {
208                 vn             = v[n][d];
209                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
210                 {
211                     vv             = lg*vn + f[n][d]*w_dt;
212
213                     /* do not scale the mean velocities u */
214                     u              = gstat[ga].u[d];
215                     va             = vv + accel[ga][d]*dt;
216                     vb             = va + (1.0-lg)*u;
217                     v[n][d]        = vb;
218                     xprime[n][d]   = x[n][d]+vb*dt;
219                 }
220                 else
221                 {
222                     v[n][d]        = 0.0;
223                     xprime[n][d]   = x[n][d];
224                 }
225             }
226         }
227     }
228     else
229     {
230         /* Plain update with Berendsen/v-rescale coupling */
231         for (n = start; n < nrend; n++)
232         {
233             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
234             {
235                 w_dt = invmass[n]*dt;
236                 if (cTC)
237                 {
238                     gt = cTC[n];
239                 }
240                 lg = tcstat[gt].lambda;
241
242                 for (d = 0; d < DIM; d++)
243                 {
244                     vn           = lg*v[n][d] + f[n][d]*w_dt;
245                     v[n][d]      = vn;
246                     xprime[n][d] = x[n][d] + vn*dt;
247                 }
248             }
249             else
250             {
251                 for (d = 0; d < DIM; d++)
252                 {
253                     v[n][d]        = 0.0;
254                     xprime[n][d]   = x[n][d];
255                 }
256             }
257         }
258     }
259 }
260
261 static void do_update_vv_vel(int start, int nrend, double dt,
262                              rvec accel[], ivec nFreeze[], real invmass[],
263                              unsigned short ptype[], unsigned short cFREEZE[],
264                              unsigned short cACC[], rvec v[], rvec f[],
265                              gmx_bool bExtended, real veta, real alpha)
266 {
267     double w_dt;
268     int    gf = 0, ga = 0;
269     int    n, d;
270     double g, mv1, mv2;
271
272     if (bExtended)
273     {
274         g        = 0.25*dt*veta*alpha;
275         mv1      = exp(-g);
276         mv2      = series_sinhx(g);
277     }
278     else
279     {
280         mv1      = 1.0;
281         mv2      = 1.0;
282     }
283     for (n = start; n < nrend; n++)
284     {
285         w_dt = invmass[n]*dt;
286         if (cFREEZE)
287         {
288             gf   = cFREEZE[n];
289         }
290         if (cACC)
291         {
292             ga   = cACC[n];
293         }
294
295         for (d = 0; d < DIM; d++)
296         {
297             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
298             {
299                 v[n][d]             = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
300             }
301             else
302             {
303                 v[n][d]        = 0.0;
304             }
305         }
306     }
307 } /* do_update_vv_vel */
308
309 static void do_update_vv_pos(int start, int nrend, double dt,
310                              ivec nFreeze[],
311                              unsigned short ptype[], unsigned short cFREEZE[],
312                              rvec x[], rvec xprime[], rvec v[],
313                              gmx_bool bExtended, real veta)
314 {
315     int    gf = 0;
316     int    n, d;
317     double g, mr1, mr2;
318
319     /* Would it make more sense if Parrinello-Rahman was put here? */
320     if (bExtended)
321     {
322         g        = 0.5*dt*veta;
323         mr1      = exp(g);
324         mr2      = series_sinhx(g);
325     }
326     else
327     {
328         mr1      = 1.0;
329         mr2      = 1.0;
330     }
331
332     for (n = start; n < nrend; n++)
333     {
334
335         if (cFREEZE)
336         {
337             gf   = cFREEZE[n];
338         }
339
340         for (d = 0; d < DIM; d++)
341         {
342             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
343             {
344                 xprime[n][d]   = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
345             }
346             else
347             {
348                 xprime[n][d]   = x[n][d];
349             }
350         }
351     }
352 } /* do_update_vv_pos */
353
354 static void do_update_visc(int start, int nrend, double dt,
355                            t_grp_tcstat *tcstat,
356                            double nh_vxi[],
357                            real invmass[],
358                            unsigned short ptype[], unsigned short cTC[],
359                            rvec x[], rvec xprime[], rvec v[],
360                            rvec f[], matrix M, matrix box, real
361                            cos_accel, real vcos,
362                            gmx_bool bNH, gmx_bool bPR)
363 {
364     double imass, w_dt;
365     int    gt = 0;
366     real   vn, vc;
367     real   lg, vxi = 0, vv;
368     real   fac, cosz;
369     rvec   vrel;
370     int    n, d;
371
372     fac = 2*M_PI/(box[ZZ][ZZ]);
373
374     if (bNH || bPR)
375     {
376         /* Update with coupling to extended ensembles, used for
377          * Nose-Hoover and Parrinello-Rahman coupling
378          */
379         for (n = start; n < nrend; n++)
380         {
381             imass = invmass[n];
382             if (cTC)
383             {
384                 gt   = cTC[n];
385             }
386             lg   = tcstat[gt].lambda;
387             cosz = cos(fac*x[n][ZZ]);
388
389             copy_rvec(v[n], vrel);
390
391             vc            = cosz*vcos;
392             vrel[XX]     -= vc;
393             if (bNH)
394             {
395                 vxi        = nh_vxi[gt];
396             }
397             for (d = 0; d < DIM; d++)
398             {
399                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
400                 {
401                     vn  = (lg*vrel[d] + dt*(imass*f[n][d] - 0.5*vxi*vrel[d]
402                                             - iprod(M[d], vrel)))/(1 + 0.5*vxi*dt);
403                     if (d == XX)
404                     {
405                         vn += vc + dt*cosz*cos_accel;
406                     }
407                     v[n][d]        = vn;
408                     xprime[n][d]   = x[n][d]+vn*dt;
409                 }
410                 else
411                 {
412                     xprime[n][d]   = x[n][d];
413                 }
414             }
415         }
416     }
417     else
418     {
419         /* Classic version of update, used with berendsen coupling */
420         for (n = start; n < nrend; n++)
421         {
422             w_dt = invmass[n]*dt;
423             if (cTC)
424             {
425                 gt   = cTC[n];
426             }
427             lg   = tcstat[gt].lambda;
428             cosz = cos(fac*x[n][ZZ]);
429
430             for (d = 0; d < DIM; d++)
431             {
432                 vn             = v[n][d];
433
434                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell))
435                 {
436                     if (d == XX)
437                     {
438                         vc           = cosz*vcos;
439                         /* Do not scale the cosine velocity profile */
440                         vv           = vc + lg*(vn - vc + f[n][d]*w_dt);
441                         /* Add the cosine accelaration profile */
442                         vv          += dt*cosz*cos_accel;
443                     }
444                     else
445                     {
446                         vv           = lg*(vn + f[n][d]*w_dt);
447                     }
448                     v[n][d]        = vv;
449                     xprime[n][d]   = x[n][d]+vv*dt;
450                 }
451                 else
452                 {
453                     v[n][d]        = 0.0;
454                     xprime[n][d]   = x[n][d];
455                 }
456             }
457         }
458     }
459 }
460
461 static gmx_stochd_t *init_stochd(t_inputrec *ir)
462 {
463     gmx_stochd_t   *sd;
464     gmx_sd_const_t *sdc;
465     int             ngtc, n;
466     real            y;
467
468     snew(sd, 1);
469
470     ngtc = ir->opts.ngtc;
471
472     if (ir->eI == eiBD)
473     {
474         snew(sd->bd_rf, ngtc);
475     }
476     else if (EI_SD(ir->eI))
477     {
478         snew(sd->sdc, ngtc);
479         snew(sd->sdsig, ngtc);
480
481         sdc = sd->sdc;
482         for (n = 0; n < ngtc; n++)
483         {
484             if (ir->opts.tau_t[n] > 0)
485             {
486                 sdc[n].gdt = ir->delta_t/ir->opts.tau_t[n];
487                 sdc[n].eph = exp(sdc[n].gdt/2);
488                 sdc[n].emh = exp(-sdc[n].gdt/2);
489                 sdc[n].em  = exp(-sdc[n].gdt);
490             }
491             else
492             {
493                 /* No friction and noise on this group */
494                 sdc[n].gdt = 0;
495                 sdc[n].eph = 1;
496                 sdc[n].emh = 1;
497                 sdc[n].em  = 1;
498             }
499             if (sdc[n].gdt >= 0.05)
500             {
501                 sdc[n].b = sdc[n].gdt*(sdc[n].eph*sdc[n].eph - 1)
502                     - 4*(sdc[n].eph - 1)*(sdc[n].eph - 1);
503                 sdc[n].c = sdc[n].gdt - 3 + 4*sdc[n].emh - sdc[n].em;
504                 sdc[n].d = 2 - sdc[n].eph - sdc[n].emh;
505             }
506             else
507             {
508                 y = sdc[n].gdt/2;
509                 /* Seventh order expansions for small y */
510                 sdc[n].b = y*y*y*y*(1/3.0+y*(1/3.0+y*(17/90.0+y*7/9.0)));
511                 sdc[n].c = y*y*y*(2/3.0+y*(-1/2.0+y*(7/30.0+y*(-1/12.0+y*31/1260.0))));
512                 sdc[n].d = y*y*(-1+y*y*(-1/12.0-y*y/360.0));
513             }
514             if (debug)
515             {
516                 fprintf(debug, "SD const tc-grp %d: b %g  c %g  d %g\n",
517                         n, sdc[n].b, sdc[n].c, sdc[n].d);
518             }
519         }
520     }
521     else if (ETC_ANDERSEN(ir->etc))
522     {
523         int        ngtc;
524         t_grpopts *opts;
525         real       reft;
526
527         opts = &ir->opts;
528         ngtc = opts->ngtc;
529
530         snew(sd->randomize_group, ngtc);
531         snew(sd->boltzfac, ngtc);
532
533         /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
534         /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
535
536         for (n = 0; n < ngtc; n++)
537         {
538             reft = std::max<real>(0, opts->ref_t[n]);
539             if ((opts->tau_t[n] > 0) && (reft > 0))  /* tau_t or ref_t = 0 means that no randomization is done */
540             {
541                 sd->randomize_group[n] = TRUE;
542                 sd->boltzfac[n]        = BOLTZ*opts->ref_t[n];
543             }
544             else
545             {
546                 sd->randomize_group[n] = FALSE;
547             }
548         }
549     }
550     return sd;
551 }
552
553 gmx_update_t init_update(t_inputrec *ir)
554 {
555     t_gmx_update *upd;
556
557     snew(upd, 1);
558
559     if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
560     {
561         upd->sd    = init_stochd(ir);
562     }
563
564     upd->xp        = NULL;
565     upd->xp_nalloc = 0;
566
567     return upd;
568 }
569
570 static void do_update_sd1(gmx_stochd_t *sd,
571                           int start, int nrend, double dt,
572                           rvec accel[], ivec nFreeze[],
573                           real invmass[], unsigned short ptype[],
574                           unsigned short cFREEZE[], unsigned short cACC[],
575                           unsigned short cTC[],
576                           rvec x[], rvec xprime[], rvec v[], rvec f[],
577                           int ngtc, real ref_t[],
578                           gmx_bool bDoConstr,
579                           gmx_bool bFirstHalfConstr,
580                           gmx_int64_t step, int seed, int* gatindex)
581 {
582     gmx_sd_const_t *sdc;
583     gmx_sd_sigma_t *sig;
584     real            kT;
585     int             gf = 0, ga = 0, gt = 0;
586     real            ism;
587     int             n, d;
588
589     sdc = sd->sdc;
590     sig = sd->sdsig;
591
592     for (n = 0; n < ngtc; n++)
593     {
594         kT = BOLTZ*ref_t[n];
595         /* The mass is encounted for later, since this differs per atom */
596         sig[n].V  = sqrt(kT*(1 - sdc[n].em*sdc[n].em));
597     }
598
599     if (!bDoConstr)
600     {
601         for (n = start; n < nrend; n++)
602         {
603             real rnd[3];
604             int  ng = gatindex ? gatindex[n] : n;
605
606             ism = sqrt(invmass[n]);
607             if (cFREEZE)
608             {
609                 gf  = cFREEZE[n];
610             }
611             if (cACC)
612             {
613                 ga  = cACC[n];
614             }
615             if (cTC)
616             {
617                 gt  = cTC[n];
618             }
619
620             gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
621
622             for (d = 0; d < DIM; d++)
623             {
624                 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
625                 {
626                     real sd_V, vn;
627
628                     sd_V         = ism*sig[gt].V*rnd[d];
629                     vn           = v[n][d] + (invmass[n]*f[n][d] + accel[ga][d])*dt;
630                     v[n][d]      = vn*sdc[gt].em + sd_V;
631                     /* Here we include half of the friction+noise
632                      * update of v into the integration of x.
633                      */
634                     xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
635                 }
636                 else
637                 {
638                     v[n][d]      = 0.0;
639                     xprime[n][d] = x[n][d];
640                 }
641             }
642         }
643     }
644     else
645     {
646         /* We do have constraints */
647         if (bFirstHalfConstr)
648         {
649             /* First update without friction and noise */
650             real im;
651
652             for (n = start; n < nrend; n++)
653             {
654                 im = invmass[n];
655
656                 if (cFREEZE)
657                 {
658                     gf  = cFREEZE[n];
659                 }
660                 if (cACC)
661                 {
662                     ga  = cACC[n];
663                 }
664
665                 for (d = 0; d < DIM; d++)
666                 {
667                     if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
668                     {
669                         v[n][d]      = v[n][d] + (im*f[n][d] + accel[ga][d])*dt;
670                         xprime[n][d] = x[n][d] +  v[n][d]*dt;
671                     }
672                     else
673                     {
674                         v[n][d]      = 0.0;
675                         xprime[n][d] = x[n][d];
676                     }
677                 }
678             }
679         }
680         else
681         {
682             /* Update friction and noise only */
683             for (n = start; n < nrend; n++)
684             {
685                 real rnd[3];
686                 int  ng = gatindex ? gatindex[n] : n;
687
688                 ism = sqrt(invmass[n]);
689                 if (cFREEZE)
690                 {
691                     gf  = cFREEZE[n];
692                 }
693                 if (cTC)
694                 {
695                     gt  = cTC[n];
696                 }
697
698                 gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
699
700                 for (d = 0; d < DIM; d++)
701                 {
702                     if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
703                     {
704                         real sd_V, vn;
705
706                         sd_V         = ism*sig[gt].V*rnd[d];
707                         vn           = v[n][d];
708                         v[n][d]      = vn*sdc[gt].em + sd_V;
709                         /* Add the friction and noise contribution only */
710                         xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
711                     }
712                 }
713             }
714         }
715     }
716 }
717
718 static void check_sd2_work_data_allocation(gmx_stochd_t *sd, int nrend)
719 {
720     if (nrend > sd->sd_V_nalloc)
721     {
722         sd->sd_V_nalloc = over_alloc_dd(nrend);
723         srenew(sd->sd_V, sd->sd_V_nalloc);
724     }
725 }
726
727 static void do_update_sd2_Tconsts(gmx_stochd_t *sd,
728                                   int           ngtc,
729                                   const real    tau_t[],
730                                   const real    ref_t[])
731 {
732     /* This is separated from the update below, because it is single threaded */
733     gmx_sd_const_t *sdc;
734     gmx_sd_sigma_t *sig;
735     int             gt;
736     real            kT;
737
738     sdc = sd->sdc;
739     sig = sd->sdsig;
740
741     for (gt = 0; gt < ngtc; gt++)
742     {
743         kT = BOLTZ*ref_t[gt];
744         /* The mass is encounted for later, since this differs per atom */
745         sig[gt].V  = sqrt(kT*(1-sdc[gt].em));
746         sig[gt].X  = sqrt(kT*sqr(tau_t[gt])*sdc[gt].c);
747         sig[gt].Yv = sqrt(kT*sdc[gt].b/sdc[gt].c);
748         sig[gt].Yx = sqrt(kT*sqr(tau_t[gt])*sdc[gt].b/(1-sdc[gt].em));
749     }
750 }
751
752 static void do_update_sd2(gmx_stochd_t *sd,
753                           gmx_bool bInitStep,
754                           int start, int nrend,
755                           rvec accel[], ivec nFreeze[],
756                           real invmass[], unsigned short ptype[],
757                           unsigned short cFREEZE[], unsigned short cACC[],
758                           unsigned short cTC[],
759                           rvec x[], rvec xprime[], rvec v[], rvec f[],
760                           rvec sd_X[],
761                           const real tau_t[],
762                           gmx_bool bFirstHalf, gmx_int64_t step, int seed,
763                           int* gatindex)
764 {
765     gmx_sd_const_t *sdc;
766     gmx_sd_sigma_t *sig;
767     /* The random part of the velocity update, generated in the first
768      * half of the update, needs to be remembered for the second half.
769      */
770     rvec  *sd_V;
771     int    gf = 0, ga = 0, gt = 0;
772     real   vn = 0, Vmh, Xmh;
773     real   ism;
774     int    n, d, ng;
775
776     sdc  = sd->sdc;
777     sig  = sd->sdsig;
778     sd_V = sd->sd_V;
779
780     for (n = start; n < nrend; n++)
781     {
782         real rnd[6], rndi[3];
783         ng  = gatindex ? gatindex[n] : n;
784         ism = sqrt(invmass[n]);
785         if (cFREEZE)
786         {
787             gf  = cFREEZE[n];
788         }
789         if (cACC)
790         {
791             ga  = cACC[n];
792         }
793         if (cTC)
794         {
795             gt  = cTC[n];
796         }
797
798         gmx_rng_cycle_6gaussian_table(step*2+(bFirstHalf ? 1 : 2), ng, seed, RND_SEED_UPDATE, rnd);
799         if (bInitStep)
800         {
801             gmx_rng_cycle_3gaussian_table(step*2, ng, seed, RND_SEED_UPDATE, rndi);
802         }
803         for (d = 0; d < DIM; d++)
804         {
805             if (bFirstHalf)
806             {
807                 vn             = v[n][d];
808             }
809             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
810             {
811                 if (bFirstHalf)
812                 {
813                     if (bInitStep)
814                     {
815                         sd_X[n][d] = ism*sig[gt].X*rndi[d];
816                     }
817                     Vmh = sd_X[n][d]*sdc[gt].d/(tau_t[gt]*sdc[gt].c)
818                         + ism*sig[gt].Yv*rnd[d*2];
819                     sd_V[n][d] = ism*sig[gt].V*rnd[d*2+1];
820
821                     v[n][d] = vn*sdc[gt].em
822                         + (invmass[n]*f[n][d] + accel[ga][d])*tau_t[gt]*(1 - sdc[gt].em)
823                         + sd_V[n][d] - sdc[gt].em*Vmh;
824
825                     xprime[n][d] = x[n][d] + v[n][d]*tau_t[gt]*(sdc[gt].eph - sdc[gt].emh);
826                 }
827                 else
828                 {
829                     /* Correct the velocities for the constraints.
830                      * This operation introduces some inaccuracy,
831                      * since the velocity is determined from differences in coordinates.
832                      */
833                     v[n][d] =
834                         (xprime[n][d] - x[n][d])/(tau_t[gt]*(sdc[gt].eph - sdc[gt].emh));
835
836                     Xmh = sd_V[n][d]*tau_t[gt]*sdc[gt].d/(sdc[gt].em-1)
837                         + ism*sig[gt].Yx*rnd[d*2];
838                     sd_X[n][d] = ism*sig[gt].X*rnd[d*2+1];
839
840                     xprime[n][d] += sd_X[n][d] - Xmh;
841
842                 }
843             }
844             else
845             {
846                 if (bFirstHalf)
847                 {
848                     v[n][d]        = 0.0;
849                     xprime[n][d]   = x[n][d];
850                 }
851             }
852         }
853     }
854 }
855
856 static void do_update_bd_Tconsts(double dt, real friction_coefficient,
857                                  int ngtc, const real ref_t[],
858                                  real *rf)
859 {
860     /* This is separated from the update below, because it is single threaded */
861     int gt;
862
863     if (friction_coefficient != 0)
864     {
865         for (gt = 0; gt < ngtc; gt++)
866         {
867             rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]/(friction_coefficient*dt));
868         }
869     }
870     else
871     {
872         for (gt = 0; gt < ngtc; gt++)
873         {
874             rf[gt] = sqrt(2.0*BOLTZ*ref_t[gt]);
875         }
876     }
877 }
878
879 static void do_update_bd(int start, int nrend, double dt,
880                          ivec nFreeze[],
881                          real invmass[], unsigned short ptype[],
882                          unsigned short cFREEZE[], unsigned short cTC[],
883                          rvec x[], rvec xprime[], rvec v[],
884                          rvec f[], real friction_coefficient,
885                          real *rf, gmx_int64_t step, int seed,
886                          int* gatindex)
887 {
888     /* note -- these appear to be full step velocities . . .  */
889     int    gf = 0, gt = 0;
890     real   vn;
891     real   invfr = 0;
892     int    n, d;
893
894     if (friction_coefficient != 0)
895     {
896         invfr = 1.0/friction_coefficient;
897     }
898
899     for (n = start; (n < nrend); n++)
900     {
901         real rnd[3];
902         int  ng  = gatindex ? gatindex[n] : n;
903
904         if (cFREEZE)
905         {
906             gf = cFREEZE[n];
907         }
908         if (cTC)
909         {
910             gt = cTC[n];
911         }
912         gmx_rng_cycle_3gaussian_table(step, ng, seed, RND_SEED_UPDATE, rnd);
913         for (d = 0; (d < DIM); d++)
914         {
915             if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
916             {
917                 if (friction_coefficient != 0)
918                 {
919                     vn = invfr*f[n][d] + rf[gt]*rnd[d];
920                 }
921                 else
922                 {
923                     /* NOTE: invmass = 2/(mass*friction_constant*dt) */
924                     vn = 0.5*invmass[n]*f[n][d]*dt
925                         + sqrt(0.5*invmass[n])*rf[gt]*rnd[d];
926                 }
927
928                 v[n][d]      = vn;
929                 xprime[n][d] = x[n][d]+vn*dt;
930             }
931             else
932             {
933                 v[n][d]      = 0.0;
934                 xprime[n][d] = x[n][d];
935             }
936         }
937     }
938 }
939
940 static void dump_it_all(FILE gmx_unused *fp, const char gmx_unused *title,
941                         int gmx_unused natoms, rvec gmx_unused x[], rvec gmx_unused xp[],
942                         rvec gmx_unused v[], rvec gmx_unused f[])
943 {
944 #ifdef DEBUG
945     if (fp)
946     {
947         fprintf(fp, "%s\n", title);
948         pr_rvecs(fp, 0, "x", x, natoms);
949         pr_rvecs(fp, 0, "xp", xp, natoms);
950         pr_rvecs(fp, 0, "v", v, natoms);
951         pr_rvecs(fp, 0, "f", f, natoms);
952     }
953 #endif
954 }
955
956 static void calc_ke_part_normal(rvec v[], t_grpopts *opts, t_mdatoms *md,
957                                 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel,
958                                 gmx_bool bSaveEkinOld)
959 {
960     int           g;
961     t_grp_tcstat *tcstat  = ekind->tcstat;
962     t_grp_acc    *grpstat = ekind->grpstat;
963     int           nthread, thread;
964
965     /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin.  Leap with AveVel is also
966        an option, but not supported now.  Additionally, if we are doing iterations.
967        bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
968        bSavEkinOld: If TRUE (in the case of iteration = bIterate is TRUE), we don't copy over the ekinh_old.
969        If FALSE, we overrwrite it.
970      */
971
972     /* group velocities are calculated in update_ekindata and
973      * accumulated in acumulate_groups.
974      * Now the partial global and groups ekin.
975      */
976     for (g = 0; (g < opts->ngtc); g++)
977     {
978
979         if (!bSaveEkinOld)
980         {
981             copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
982         }
983         if (bEkinAveVel)
984         {
985             clear_mat(tcstat[g].ekinf);
986         }
987         else
988         {
989             clear_mat(tcstat[g].ekinh);
990         }
991         if (bEkinAveVel)
992         {
993             tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
994         }
995     }
996     ekind->dekindl_old = ekind->dekindl;
997
998     nthread = gmx_omp_nthreads_get(emntUpdate);
999
1000 #pragma omp parallel for num_threads(nthread) schedule(static)
1001     for (thread = 0; thread < nthread; thread++)
1002     {
1003         int     start_t, end_t, n;
1004         int     ga, gt;
1005         rvec    v_corrt;
1006         real    hm;
1007         int     d, m;
1008         matrix *ekin_sum;
1009         real   *dekindl_sum;
1010
1011         start_t = ((thread+0)*md->homenr)/nthread;
1012         end_t   = ((thread+1)*md->homenr)/nthread;
1013
1014         ekin_sum    = ekind->ekin_work[thread];
1015         dekindl_sum = ekind->dekindl_work[thread];
1016
1017         for (gt = 0; gt < opts->ngtc; gt++)
1018         {
1019             clear_mat(ekin_sum[gt]);
1020         }
1021         *dekindl_sum = 0.0;
1022
1023         ga = 0;
1024         gt = 0;
1025         for (n = start_t; n < end_t; n++)
1026         {
1027             if (md->cACC)
1028             {
1029                 ga = md->cACC[n];
1030             }
1031             if (md->cTC)
1032             {
1033                 gt = md->cTC[n];
1034             }
1035             hm   = 0.5*md->massT[n];
1036
1037             for (d = 0; (d < DIM); d++)
1038             {
1039                 v_corrt[d]  = v[n][d]  - grpstat[ga].u[d];
1040             }
1041             for (d = 0; (d < DIM); d++)
1042             {
1043                 for (m = 0; (m < DIM); m++)
1044                 {
1045                     /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
1046                     ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1047                 }
1048             }
1049             if (md->nMassPerturbed && md->bPerturbed[n])
1050             {
1051                 *dekindl_sum +=
1052                     0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1053             }
1054         }
1055     }
1056
1057     ekind->dekindl = 0;
1058     for (thread = 0; thread < nthread; thread++)
1059     {
1060         for (g = 0; g < opts->ngtc; g++)
1061         {
1062             if (bEkinAveVel)
1063             {
1064                 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1065                       tcstat[g].ekinf);
1066             }
1067             else
1068             {
1069                 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1070                       tcstat[g].ekinh);
1071             }
1072         }
1073
1074         ekind->dekindl += *ekind->dekindl_work[thread];
1075     }
1076
1077     inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1078 }
1079
1080 static void calc_ke_part_visc(matrix box, rvec x[], rvec v[],
1081                               t_grpopts *opts, t_mdatoms *md,
1082                               gmx_ekindata_t *ekind,
1083                               t_nrnb *nrnb, gmx_bool bEkinAveVel)
1084 {
1085     int           start = 0, homenr = md->homenr;
1086     int           g, d, n, m, gt = 0;
1087     rvec          v_corrt;
1088     real          hm;
1089     t_grp_tcstat *tcstat = ekind->tcstat;
1090     t_cos_acc    *cosacc = &(ekind->cosacc);
1091     real          dekindl;
1092     real          fac, cosz;
1093     double        mvcos;
1094
1095     for (g = 0; g < opts->ngtc; g++)
1096     {
1097         copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1098         clear_mat(ekind->tcstat[g].ekinh);
1099     }
1100     ekind->dekindl_old = ekind->dekindl;
1101
1102     fac     = 2*M_PI/box[ZZ][ZZ];
1103     mvcos   = 0;
1104     dekindl = 0;
1105     for (n = start; n < start+homenr; n++)
1106     {
1107         if (md->cTC)
1108         {
1109             gt = md->cTC[n];
1110         }
1111         hm   = 0.5*md->massT[n];
1112
1113         /* Note that the times of x and v differ by half a step */
1114         /* MRS -- would have to be changed for VV */
1115         cosz         = cos(fac*x[n][ZZ]);
1116         /* Calculate the amplitude of the new velocity profile */
1117         mvcos       += 2*cosz*md->massT[n]*v[n][XX];
1118
1119         copy_rvec(v[n], v_corrt);
1120         /* Subtract the profile for the kinetic energy */
1121         v_corrt[XX] -= cosz*cosacc->vcos;
1122         for (d = 0; (d < DIM); d++)
1123         {
1124             for (m = 0; (m < DIM); m++)
1125             {
1126                 /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
1127                 if (bEkinAveVel)
1128                 {
1129                     tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1130                 }
1131                 else
1132                 {
1133                     tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1134                 }
1135             }
1136         }
1137         if (md->nPerturbed && md->bPerturbed[n])
1138         {
1139             /* The minus sign here might be confusing.
1140              * The kinetic contribution from dH/dl doesn't come from
1141              * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1142              * where p are the momenta. The difference is only a minus sign.
1143              */
1144             dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1145         }
1146     }
1147     ekind->dekindl = dekindl;
1148     cosacc->mvcos  = mvcos;
1149
1150     inc_nrnb(nrnb, eNR_EKIN, homenr);
1151 }
1152
1153 void calc_ke_part(t_state *state, t_grpopts *opts, t_mdatoms *md,
1154                   gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel, gmx_bool bSaveEkinOld)
1155 {
1156     if (ekind->cosacc.cos_accel == 0)
1157     {
1158         calc_ke_part_normal(state->v, opts, md, ekind, nrnb, bEkinAveVel, bSaveEkinOld);
1159     }
1160     else
1161     {
1162         calc_ke_part_visc(state->box, state->x, state->v, opts, md, ekind, nrnb, bEkinAveVel);
1163     }
1164 }
1165
1166 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1167 {
1168     ekinstate->ekin_n = ir->opts.ngtc;
1169     snew(ekinstate->ekinh, ekinstate->ekin_n);
1170     snew(ekinstate->ekinf, ekinstate->ekin_n);
1171     snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1172     snew(ekinstate->ekinscalef_nhc, ekinstate->ekin_n);
1173     snew(ekinstate->ekinscaleh_nhc, ekinstate->ekin_n);
1174     snew(ekinstate->vscale_nhc, ekinstate->ekin_n);
1175     ekinstate->dekindl = 0;
1176     ekinstate->mvcos   = 0;
1177 }
1178
1179 void update_ekinstate(ekinstate_t *ekinstate, gmx_ekindata_t *ekind)
1180 {
1181     int i;
1182
1183     for (i = 0; i < ekinstate->ekin_n; i++)
1184     {
1185         copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1186         copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1187         copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1188         ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1189         ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1190         ekinstate->vscale_nhc[i]     = ekind->tcstat[i].vscale_nhc;
1191     }
1192
1193     copy_mat(ekind->ekin, ekinstate->ekin_total);
1194     ekinstate->dekindl = ekind->dekindl;
1195     ekinstate->mvcos   = ekind->cosacc.mvcos;
1196
1197 }
1198
1199 void restore_ekinstate_from_state(t_commrec *cr,
1200                                   gmx_ekindata_t *ekind, ekinstate_t *ekinstate)
1201 {
1202     int i, n;
1203
1204     if (MASTER(cr))
1205     {
1206         for (i = 0; i < ekinstate->ekin_n; i++)
1207         {
1208             copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1209             copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1210             copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1211             ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1212             ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1213             ekind->tcstat[i].vscale_nhc     = ekinstate->vscale_nhc[i];
1214         }
1215
1216         copy_mat(ekinstate->ekin_total, ekind->ekin);
1217
1218         ekind->dekindl      = ekinstate->dekindl;
1219         ekind->cosacc.mvcos = ekinstate->mvcos;
1220         n                   = ekinstate->ekin_n;
1221     }
1222
1223     if (PAR(cr))
1224     {
1225         gmx_bcast(sizeof(n), &n, cr);
1226         for (i = 0; i < n; i++)
1227         {
1228             gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1229                       ekind->tcstat[i].ekinh[0], cr);
1230             gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1231                       ekind->tcstat[i].ekinf[0], cr);
1232             gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1233                       ekind->tcstat[i].ekinh_old[0], cr);
1234
1235             gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1236                       &(ekind->tcstat[i].ekinscalef_nhc), cr);
1237             gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1238                       &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1239             gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1240                       &(ekind->tcstat[i].vscale_nhc), cr);
1241         }
1242         gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1243                   ekind->ekin[0], cr);
1244
1245         gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1246         gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1247     }
1248 }
1249
1250 void set_deform_reference_box(gmx_update_t upd, gmx_int64_t step, matrix box)
1251 {
1252     upd->deformref_step = step;
1253     copy_mat(box, upd->deformref_box);
1254 }
1255
1256 static void deform(gmx_update_t upd,
1257                    int start, int homenr, rvec x[], matrix box, matrix *scale_tot,
1258                    const t_inputrec *ir, gmx_int64_t step)
1259 {
1260     matrix bnew, invbox, mu;
1261     real   elapsed_time;
1262     int    i, j;
1263
1264     elapsed_time = (step + 1 - upd->deformref_step)*ir->delta_t;
1265     copy_mat(box, bnew);
1266     for (i = 0; i < DIM; i++)
1267     {
1268         for (j = 0; j < DIM; j++)
1269         {
1270             if (ir->deform[i][j] != 0)
1271             {
1272                 bnew[i][j] =
1273                     upd->deformref_box[i][j] + elapsed_time*ir->deform[i][j];
1274             }
1275         }
1276     }
1277     /* We correct the off-diagonal elements,
1278      * which can grow indefinitely during shearing,
1279      * so the shifts do not get messed up.
1280      */
1281     for (i = 1; i < DIM; i++)
1282     {
1283         for (j = i-1; j >= 0; j--)
1284         {
1285             while (bnew[i][j] - box[i][j] > 0.5*bnew[j][j])
1286             {
1287                 rvec_dec(bnew[i], bnew[j]);
1288             }
1289             while (bnew[i][j] - box[i][j] < -0.5*bnew[j][j])
1290             {
1291                 rvec_inc(bnew[i], bnew[j]);
1292             }
1293         }
1294     }
1295     m_inv_ur0(box, invbox);
1296     copy_mat(bnew, box);
1297     mmul_ur0(box, invbox, mu);
1298
1299     for (i = start; i < start+homenr; i++)
1300     {
1301         x[i][XX] = mu[XX][XX]*x[i][XX]+mu[YY][XX]*x[i][YY]+mu[ZZ][XX]*x[i][ZZ];
1302         x[i][YY] = mu[YY][YY]*x[i][YY]+mu[ZZ][YY]*x[i][ZZ];
1303         x[i][ZZ] = mu[ZZ][ZZ]*x[i][ZZ];
1304     }
1305     if (scale_tot != NULL)
1306     {
1307         /* The transposes of the scaling matrices are stored,
1308          * so we need to do matrix multiplication in the inverse order.
1309          */
1310         mmul_ur0(*scale_tot, mu, *scale_tot);
1311     }
1312 }
1313
1314 void update_tcouple(gmx_int64_t       step,
1315                     t_inputrec       *inputrec,
1316                     t_state          *state,
1317                     gmx_ekindata_t   *ekind,
1318                     t_extmass        *MassQ,
1319                     t_mdatoms        *md)
1320
1321 {
1322     gmx_bool   bTCouple = FALSE;
1323     real       dttc;
1324     int        i, offset;
1325
1326     /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1327     if (inputrec->etc != etcNO &&
1328         !(IR_NVT_TROTTER(inputrec) || IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec)))
1329     {
1330         /* We should only couple after a step where energies were determined (for leapfrog versions)
1331            or the step energies are determined, for velocity verlet versions */
1332
1333         if (EI_VV(inputrec->eI))
1334         {
1335             offset = 0;
1336         }
1337         else
1338         {
1339             offset = 1;
1340         }
1341         bTCouple = (inputrec->nsttcouple == 1 ||
1342                     do_per_step(step+inputrec->nsttcouple-offset,
1343                                 inputrec->nsttcouple));
1344     }
1345
1346     if (bTCouple)
1347     {
1348         dttc = inputrec->nsttcouple*inputrec->delta_t;
1349
1350         switch (inputrec->etc)
1351         {
1352             case etcNO:
1353                 break;
1354             case etcBERENDSEN:
1355                 berendsen_tcoupl(inputrec, ekind, dttc);
1356                 break;
1357             case etcNOSEHOOVER:
1358                 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1359                                   state->nosehoover_xi, state->nosehoover_vxi, MassQ);
1360                 break;
1361             case etcVRESCALE:
1362                 vrescale_tcoupl(inputrec, step, ekind, dttc,
1363                                 state->therm_integral);
1364                 break;
1365         }
1366         /* rescale in place here */
1367         if (EI_VV(inputrec->eI))
1368         {
1369             rescale_velocities(ekind, md, 0, md->homenr, state->v);
1370         }
1371     }
1372     else
1373     {
1374         /* Set the T scaling lambda to 1 to have no scaling */
1375         for (i = 0; (i < inputrec->opts.ngtc); i++)
1376         {
1377             ekind->tcstat[i].lambda = 1.0;
1378         }
1379     }
1380 }
1381
1382 void update_pcouple(FILE             *fplog,
1383                     gmx_int64_t       step,
1384                     t_inputrec       *inputrec,
1385                     t_state          *state,
1386                     matrix            pcoupl_mu,
1387                     matrix            M,
1388                     gmx_bool          bInitStep)
1389 {
1390     gmx_bool   bPCouple = FALSE;
1391     real       dtpc     = 0;
1392     int        i;
1393
1394     /* if using Trotter pressure, we do this in coupling.c, so we leave it false. */
1395     if (inputrec->epc != epcNO && (!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))))
1396     {
1397         /* We should only couple after a step where energies were determined */
1398         bPCouple = (inputrec->nstpcouple == 1 ||
1399                     do_per_step(step+inputrec->nstpcouple-1,
1400                                 inputrec->nstpcouple));
1401     }
1402
1403     clear_mat(pcoupl_mu);
1404     for (i = 0; i < DIM; i++)
1405     {
1406         pcoupl_mu[i][i] = 1.0;
1407     }
1408
1409     clear_mat(M);
1410
1411     if (bPCouple)
1412     {
1413         dtpc = inputrec->nstpcouple*inputrec->delta_t;
1414
1415         switch (inputrec->epc)
1416         {
1417             /* We can always pcoupl, even if we did not sum the energies
1418              * the previous step, since state->pres_prev is only updated
1419              * when the energies have been summed.
1420              */
1421             case (epcNO):
1422                 break;
1423             case (epcBERENDSEN):
1424                 if (!bInitStep)
1425                 {
1426                     berendsen_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev, state->box,
1427                                      pcoupl_mu);
1428                 }
1429                 break;
1430             case (epcPARRINELLORAHMAN):
1431                 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1432                                         state->box, state->box_rel, state->boxv,
1433                                         M, pcoupl_mu, bInitStep);
1434                 break;
1435             default:
1436                 break;
1437         }
1438     }
1439 }
1440
1441 static rvec *get_xprime(const t_state *state, gmx_update_t upd)
1442 {
1443     if (state->nalloc > upd->xp_nalloc)
1444     {
1445         upd->xp_nalloc = state->nalloc;
1446         srenew(upd->xp, upd->xp_nalloc);
1447     }
1448
1449     return upd->xp;
1450 }
1451
1452 static void combine_forces(gmx_update_t upd,
1453                            int nstcalclr,
1454                            gmx_constr_t constr,
1455                            t_inputrec *ir, t_mdatoms *md, t_idef *idef,
1456                            t_commrec *cr,
1457                            gmx_int64_t step,
1458                            t_state *state, gmx_bool bMolPBC,
1459                            int start, int nrend,
1460                            rvec f[], rvec f_lr[],
1461                            tensor *vir_lr_constr,
1462                            t_nrnb *nrnb)
1463 {
1464     int  i, d;
1465
1466     /* f contains the short-range forces + the long range forces
1467      * which are stored separately in f_lr.
1468      */
1469
1470     if (constr != NULL && vir_lr_constr != NULL &&
1471         !(ir->eConstrAlg == econtSHAKE && ir->epc == epcNO))
1472     {
1473         /* We need to constrain the LR forces separately,
1474          * because due to the different pre-factor for the SR and LR
1475          * forces in the update algorithm, we have to correct
1476          * the constraint virial for the nstcalclr-1 extra f_lr.
1477          * Constrain only the additional LR part of the force.
1478          */
1479         /* MRS -- need to make sure this works with trotter integration -- the constraint calls may not be right.*/
1480         rvec *xp;
1481         real  fac;
1482         int   gf = 0;
1483
1484         xp  = get_xprime(state, upd);
1485
1486         fac = (nstcalclr - 1)*ir->delta_t*ir->delta_t;
1487
1488         for (i = 0; i < md->homenr; i++)
1489         {
1490             if (md->cFREEZE != NULL)
1491             {
1492                 gf = md->cFREEZE[i];
1493             }
1494             for (d = 0; d < DIM; d++)
1495             {
1496                 if ((md->ptype[i] != eptVSite) &&
1497                     (md->ptype[i] != eptShell) &&
1498                     !ir->opts.nFreeze[gf][d])
1499                 {
1500                     xp[i][d] = state->x[i][d] + fac*f_lr[i][d]*md->invmass[i];
1501                 }
1502                 else
1503                 {
1504                     xp[i][d] = state->x[i][d];
1505                 }
1506             }
1507         }
1508         constrain(NULL, FALSE, FALSE, constr, idef, ir, NULL, cr, step, 0, 1.0, md,
1509                   state->x, xp, xp, bMolPBC, state->box, state->lambda[efptBONDED], NULL,
1510                   NULL, vir_lr_constr, nrnb, econqCoord, ir->epc == epcMTTK, state->veta, state->veta);
1511     }
1512
1513     /* Add nstcalclr-1 times the LR force to the sum of both forces
1514      * and store the result in forces_lr.
1515      */
1516     for (i = start; i < nrend; i++)
1517     {
1518         for (d = 0; d < DIM; d++)
1519         {
1520             f_lr[i][d] = f[i][d] + (nstcalclr - 1)*f_lr[i][d];
1521         }
1522     }
1523 }
1524
1525 void update_constraints(FILE             *fplog,
1526                         gmx_int64_t       step,
1527                         real             *dvdlambda, /* the contribution to be added to the bonded interactions */
1528                         t_inputrec       *inputrec,  /* input record and box stuff      */
1529                         gmx_ekindata_t   *ekind,
1530                         t_mdatoms        *md,
1531                         t_state          *state,
1532                         gmx_bool          bMolPBC,
1533                         t_graph          *graph,
1534                         rvec              force[],   /* forces on home particles */
1535                         t_idef           *idef,
1536                         tensor            vir_part,
1537                         t_commrec        *cr,
1538                         t_nrnb           *nrnb,
1539                         gmx_wallcycle_t   wcycle,
1540                         gmx_update_t      upd,
1541                         gmx_constr_t      constr,
1542                         gmx_bool          bFirstHalf,
1543                         gmx_bool          bCalcVir,
1544                         real              vetanew)
1545 {
1546     gmx_bool             bLastStep, bLog = FALSE, bEner = FALSE, bDoConstr = FALSE;
1547     double               dt;
1548     int                  start, homenr, nrend, i, m;
1549     tensor               vir_con;
1550     rvec                *xprime = NULL;
1551     int                  nth, th;
1552
1553     if (constr)
1554     {
1555         bDoConstr = TRUE;
1556     }
1557     if (bFirstHalf && !EI_VV(inputrec->eI))
1558     {
1559         bDoConstr = FALSE;
1560     }
1561
1562     /* for now, SD update is here -- though it really seems like it
1563        should be reformulated as a velocity verlet method, since it has two parts */
1564
1565     start  = 0;
1566     homenr = md->homenr;
1567     nrend  = start+homenr;
1568
1569     dt   = inputrec->delta_t;
1570
1571     /*
1572      *  Steps (7C, 8C)
1573      *  APPLY CONSTRAINTS:
1574      *  BLOCK SHAKE
1575
1576      * When doing PR pressure coupling we have to constrain the
1577      * bonds in each iteration. If we are only using Nose-Hoover tcoupling
1578      * it is enough to do this once though, since the relative velocities
1579      * after this will be normal to the bond vector
1580      */
1581
1582     if (bDoConstr)
1583     {
1584         /* clear out constraints before applying */
1585         clear_mat(vir_part);
1586
1587         xprime = get_xprime(state, upd);
1588
1589         bLastStep = (step == inputrec->init_step+inputrec->nsteps);
1590         bLog      = (do_per_step(step, inputrec->nstlog) || bLastStep || (step < 0));
1591         bEner     = (do_per_step(step, inputrec->nstenergy) || bLastStep);
1592         /* Constrain the coordinates xprime */
1593         wallcycle_start(wcycle, ewcCONSTR);
1594         if (EI_VV(inputrec->eI) && bFirstHalf)
1595         {
1596             constrain(NULL, bLog, bEner, constr, idef,
1597                       inputrec, ekind, cr, step, 1, 1.0, md,
1598                       state->x, state->v, state->v,
1599                       bMolPBC, state->box,
1600                       state->lambda[efptBONDED], dvdlambda,
1601                       NULL, bCalcVir ? &vir_con : NULL, nrnb, econqVeloc,
1602                       inputrec->epc == epcMTTK, state->veta, vetanew);
1603         }
1604         else
1605         {
1606             constrain(NULL, bLog, bEner, constr, idef,
1607                       inputrec, ekind, cr, step, 1, 1.0, md,
1608                       state->x, xprime, NULL,
1609                       bMolPBC, state->box,
1610                       state->lambda[efptBONDED], dvdlambda,
1611                       state->v, bCalcVir ? &vir_con : NULL, nrnb, econqCoord,
1612                       inputrec->epc == epcMTTK, state->veta, state->veta);
1613         }
1614         wallcycle_stop(wcycle, ewcCONSTR);
1615
1616         where();
1617
1618         dump_it_all(fplog, "After Shake",
1619                     state->natoms, state->x, xprime, state->v, force);
1620
1621         if (bCalcVir)
1622         {
1623             if (inputrec->eI == eiSD2)
1624             {
1625                 /* A correction factor eph is needed for the SD constraint force */
1626                 /* Here we can, unfortunately, not have proper corrections
1627                  * for different friction constants, so we use the first one.
1628                  */
1629                 for (i = 0; i < DIM; i++)
1630                 {
1631                     for (m = 0; m < DIM; m++)
1632                     {
1633                         vir_part[i][m] += upd->sd->sdc[0].eph*vir_con[i][m];
1634                     }
1635                 }
1636             }
1637             else
1638             {
1639                 m_add(vir_part, vir_con, vir_part);
1640             }
1641             if (debug)
1642             {
1643                 pr_rvecs(debug, 0, "constraint virial", vir_part, DIM);
1644             }
1645         }
1646     }
1647
1648     where();
1649
1650     if (inputrec->eI == eiSD1 && bDoConstr && !bFirstHalf)
1651     {
1652         wallcycle_start(wcycle, ewcUPDATE);
1653         xprime = get_xprime(state, upd);
1654
1655         nth = gmx_omp_nthreads_get(emntUpdate);
1656
1657 #pragma omp parallel for num_threads(nth) schedule(static)
1658
1659         for (th = 0; th < nth; th++)
1660         {
1661             int start_th, end_th;
1662
1663             start_th = start + ((nrend-start)* th   )/nth;
1664             end_th   = start + ((nrend-start)*(th+1))/nth;
1665
1666             /* The second part of the SD integration */
1667             do_update_sd1(upd->sd,
1668                           start_th, end_th, dt,
1669                           inputrec->opts.acc, inputrec->opts.nFreeze,
1670                           md->invmass, md->ptype,
1671                           md->cFREEZE, md->cACC, md->cTC,
1672                           state->x, xprime, state->v, force,
1673                           inputrec->opts.ngtc, inputrec->opts.ref_t,
1674                           bDoConstr, FALSE,
1675                           step, inputrec->ld_seed,
1676                           DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1677         }
1678         inc_nrnb(nrnb, eNR_UPDATE, homenr);
1679         wallcycle_stop(wcycle, ewcUPDATE);
1680
1681         if (bDoConstr)
1682         {
1683             /* Constrain the coordinates xprime for half a time step */
1684             wallcycle_start(wcycle, ewcCONSTR);
1685
1686             constrain(NULL, bLog, bEner, constr, idef,
1687                       inputrec, NULL, cr, step, 1, 0.5, md,
1688                       state->x, xprime, NULL,
1689                       bMolPBC, state->box,
1690                       state->lambda[efptBONDED], dvdlambda,
1691                       state->v, NULL, nrnb, econqCoord, FALSE, 0, 0);
1692
1693             wallcycle_stop(wcycle, ewcCONSTR);
1694         }
1695     }
1696
1697     if ((inputrec->eI == eiSD2) && !(bFirstHalf))
1698     {
1699         wallcycle_start(wcycle, ewcUPDATE);
1700         xprime = get_xprime(state, upd);
1701
1702         nth = gmx_omp_nthreads_get(emntUpdate);
1703
1704 #pragma omp parallel for num_threads(nth) schedule(static)
1705         for (th = 0; th < nth; th++)
1706         {
1707             int start_th, end_th;
1708
1709             start_th = start + ((nrend-start)* th   )/nth;
1710             end_th   = start + ((nrend-start)*(th+1))/nth;
1711
1712             /* The second part of the SD integration */
1713             do_update_sd2(upd->sd,
1714                           FALSE, start_th, end_th,
1715                           inputrec->opts.acc, inputrec->opts.nFreeze,
1716                           md->invmass, md->ptype,
1717                           md->cFREEZE, md->cACC, md->cTC,
1718                           state->x, xprime, state->v, force, state->sd_X,
1719                           inputrec->opts.tau_t,
1720                           FALSE, step, inputrec->ld_seed,
1721                           DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
1722         }
1723         inc_nrnb(nrnb, eNR_UPDATE, homenr);
1724         wallcycle_stop(wcycle, ewcUPDATE);
1725
1726         if (bDoConstr)
1727         {
1728             /* Constrain the coordinates xprime */
1729             wallcycle_start(wcycle, ewcCONSTR);
1730             constrain(NULL, bLog, bEner, constr, idef,
1731                       inputrec, NULL, cr, step, 1, 1.0, md,
1732                       state->x, xprime, NULL,
1733                       bMolPBC, state->box,
1734                       state->lambda[efptBONDED], dvdlambda,
1735                       NULL, NULL, nrnb, econqCoord, FALSE, 0, 0);
1736             wallcycle_stop(wcycle, ewcCONSTR);
1737         }
1738     }
1739
1740
1741     /* We must always unshift after updating coordinates; if we did not shake
1742        x was shifted in do_force */
1743
1744     if (!(bFirstHalf)) /* in the first half of vv, no shift. */
1745     {
1746         if (graph && (graph->nnodes > 0))
1747         {
1748             unshift_x(graph, state->box, state->x, upd->xp);
1749             if (TRICLINIC(state->box))
1750             {
1751                 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1752             }
1753             else
1754             {
1755                 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1756             }
1757         }
1758         else
1759         {
1760             // cppcheck-suppress unreadVariable
1761             nth = gmx_omp_nthreads_get(emntUpdate);
1762
1763 #pragma omp parallel for num_threads(nth) schedule(static)
1764             for (i = start; i < nrend; i++)
1765             {
1766                 copy_rvec(upd->xp[i], state->x[i]);
1767             }
1768         }
1769
1770         dump_it_all(fplog, "After unshift",
1771                     state->natoms, state->x, upd->xp, state->v, force);
1772     }
1773 /* ############# END the update of velocities and positions ######### */
1774 }
1775
1776 void update_box(FILE             *fplog,
1777                 gmx_int64_t       step,
1778                 t_inputrec       *inputrec,  /* input record and box stuff      */
1779                 t_mdatoms        *md,
1780                 t_state          *state,
1781                 rvec              force[],   /* forces on home particles */
1782                 matrix           *scale_tot,
1783                 matrix            pcoupl_mu,
1784                 t_nrnb           *nrnb,
1785                 gmx_update_t      upd)
1786 {
1787     double               dt;
1788     int                  start, homenr, i, n, m;
1789
1790     start  = 0;
1791     homenr = md->homenr;
1792
1793     dt = inputrec->delta_t;
1794
1795     where();
1796
1797     /* now update boxes */
1798     switch (inputrec->epc)
1799     {
1800         case (epcNO):
1801             break;
1802         case (epcBERENDSEN):
1803             berendsen_pscale(inputrec, pcoupl_mu, state->box, state->box_rel,
1804                              start, homenr, state->x, md->cFREEZE, nrnb);
1805             break;
1806         case (epcPARRINELLORAHMAN):
1807             /* The box velocities were updated in do_pr_pcoupl in the update
1808              * iteration, but we dont change the box vectors until we get here
1809              * since we need to be able to shift/unshift above.
1810              */
1811             for (i = 0; i < DIM; i++)
1812             {
1813                 for (m = 0; m <= i; m++)
1814                 {
1815                     state->box[i][m] += dt*state->boxv[i][m];
1816                 }
1817             }
1818             preserve_box_shape(inputrec, state->box_rel, state->box);
1819
1820             /* Scale the coordinates */
1821             for (n = start; (n < start+homenr); n++)
1822             {
1823                 tmvmul_ur0(pcoupl_mu, state->x[n], state->x[n]);
1824             }
1825             break;
1826         case (epcMTTK):
1827             switch (inputrec->epct)
1828             {
1829                 case (epctISOTROPIC):
1830                     /* DIM * eta = ln V.  so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1831                        ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1832                        Side length scales as exp(veta*dt) */
1833
1834                     msmul(state->box, exp(state->veta*dt), state->box);
1835
1836                     /* Relate veta to boxv.  veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1837                        o               If we assume isotropic scaling, and box length scaling
1838                        factor L, then V = L^DIM (det(M)).  So dV/dt = DIM
1839                        L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt.  The
1840                        determinant of B is L^DIM det(M), and the determinant
1841                        of dB/dt is (dL/dT)^DIM det (M).  veta will be
1842                        (det(dB/dT)/det(B))^(1/3).  Then since M =
1843                        B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1844
1845                     msmul(state->box, state->veta, state->boxv);
1846                     break;
1847                 default:
1848                     break;
1849             }
1850             break;
1851         default:
1852             break;
1853     }
1854
1855     if ((!(IR_NPT_TROTTER(inputrec) || IR_NPH_TROTTER(inputrec))) && scale_tot)
1856     {
1857         /* The transposes of the scaling matrices are stored,
1858          * therefore we need to reverse the order in the multiplication.
1859          */
1860         mmul_ur0(*scale_tot, pcoupl_mu, *scale_tot);
1861     }
1862
1863     if (DEFORM(*inputrec))
1864     {
1865         deform(upd, start, homenr, state->x, state->box, scale_tot, inputrec, step);
1866     }
1867     where();
1868     dump_it_all(fplog, "After update",
1869                 state->natoms, state->x, upd->xp, state->v, force);
1870 }
1871
1872 void update_coords(FILE             *fplog,
1873                    gmx_int64_t       step,
1874                    t_inputrec       *inputrec,  /* input record and box stuff   */
1875                    t_mdatoms        *md,
1876                    t_state          *state,
1877                    gmx_bool          bMolPBC,
1878                    rvec             *f,    /* forces on home particles */
1879                    gmx_bool          bDoLR,
1880                    rvec             *f_lr,
1881                    tensor           *vir_lr_constr,
1882                    t_fcdata         *fcd,
1883                    gmx_ekindata_t   *ekind,
1884                    matrix            M,
1885                    gmx_update_t      upd,
1886                    gmx_bool          bInitStep,
1887                    int               UpdatePart,
1888                    t_commrec        *cr, /* these shouldn't be here -- need to think about it */
1889                    t_nrnb           *nrnb,
1890                    gmx_constr_t      constr,
1891                    t_idef           *idef)
1892 {
1893     gmx_bool          bNH, bPR, bDoConstr = FALSE;
1894     double            dt, alpha;
1895     rvec             *force;
1896     int               start, homenr, nrend;
1897     rvec             *xprime;
1898     int               nth, th;
1899
1900     bDoConstr = (NULL != constr);
1901
1902     /* Running the velocity half does nothing except for velocity verlet */
1903     if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1904         !EI_VV(inputrec->eI))
1905     {
1906         gmx_incons("update_coords called for velocity without VV integrator");
1907     }
1908
1909     start  = 0;
1910     homenr = md->homenr;
1911     nrend  = start+homenr;
1912
1913     xprime = get_xprime(state, upd);
1914
1915     dt   = inputrec->delta_t;
1916
1917     /* We need to update the NMR restraint history when time averaging is used */
1918     if (state->flags & (1<<estDISRE_RM3TAV))
1919     {
1920         update_disres_history(fcd, &state->hist);
1921     }
1922     if (state->flags & (1<<estORIRE_DTAV))
1923     {
1924         update_orires_history(fcd, &state->hist);
1925     }
1926
1927
1928     bNH = inputrec->etc == etcNOSEHOOVER;
1929     bPR = ((inputrec->epc == epcPARRINELLORAHMAN) || (inputrec->epc == epcMTTK));
1930
1931     if (bDoLR && inputrec->nstcalclr > 1 && !EI_VV(inputrec->eI))  /* get this working with VV? */
1932     {
1933         /* Store the total force + nstcalclr-1 times the LR force
1934          * in forces_lr, so it can be used in a normal update algorithm
1935          * to produce twin time stepping.
1936          */
1937         /* is this correct in the new construction? MRS */
1938         combine_forces(upd,
1939                        inputrec->nstcalclr, constr, inputrec, md, idef, cr,
1940                        step, state, bMolPBC,
1941                        start, nrend, f, f_lr, vir_lr_constr, nrnb);
1942         force = f_lr;
1943     }
1944     else
1945     {
1946         force = f;
1947     }
1948
1949     /* ############# START The update of velocities and positions ######### */
1950     where();
1951     dump_it_all(fplog, "Before update",
1952                 state->natoms, state->x, xprime, state->v, force);
1953
1954     if (inputrec->eI == eiSD2)
1955     {
1956         check_sd2_work_data_allocation(upd->sd, nrend);
1957
1958         do_update_sd2_Tconsts(upd->sd,
1959                               inputrec->opts.ngtc,
1960                               inputrec->opts.tau_t,
1961                               inputrec->opts.ref_t);
1962     }
1963     if (inputrec->eI == eiBD)
1964     {
1965         do_update_bd_Tconsts(dt, inputrec->bd_fric,
1966                              inputrec->opts.ngtc, inputrec->opts.ref_t,
1967                              upd->sd->bd_rf);
1968     }
1969
1970     nth = gmx_omp_nthreads_get(emntUpdate);
1971
1972 #pragma omp parallel for num_threads(nth) schedule(static) private(alpha)
1973     for (th = 0; th < nth; th++)
1974     {
1975         int start_th, end_th;
1976
1977         start_th = start + ((nrend-start)* th   )/nth;
1978         end_th   = start + ((nrend-start)*(th+1))/nth;
1979
1980         switch (inputrec->eI)
1981         {
1982             case (eiMD):
1983                 if (ekind->cosacc.cos_accel == 0)
1984                 {
1985                     do_update_md(start_th, end_th, dt,
1986                                  ekind->tcstat, state->nosehoover_vxi,
1987                                  ekind->bNEMD, ekind->grpstat, inputrec->opts.acc,
1988                                  inputrec->opts.nFreeze,
1989                                  md->invmass, md->ptype,
1990                                  md->cFREEZE, md->cACC, md->cTC,
1991                                  state->x, xprime, state->v, force, M,
1992                                  bNH, bPR);
1993                 }
1994                 else
1995                 {
1996                     do_update_visc(start_th, end_th, dt,
1997                                    ekind->tcstat, state->nosehoover_vxi,
1998                                    md->invmass, md->ptype,
1999                                    md->cTC, state->x, xprime, state->v, force, M,
2000                                    state->box,
2001                                    ekind->cosacc.cos_accel,
2002                                    ekind->cosacc.vcos,
2003                                    bNH, bPR);
2004                 }
2005                 break;
2006             case (eiSD1):
2007                 /* With constraints, the SD1 update is done in 2 parts */
2008                 do_update_sd1(upd->sd,
2009                               start_th, end_th, dt,
2010                               inputrec->opts.acc, inputrec->opts.nFreeze,
2011                               md->invmass, md->ptype,
2012                               md->cFREEZE, md->cACC, md->cTC,
2013                               state->x, xprime, state->v, force,
2014                               inputrec->opts.ngtc, inputrec->opts.ref_t,
2015                               bDoConstr, TRUE,
2016                               step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2017                 break;
2018             case (eiSD2):
2019                 /* The SD2 update is always done in 2 parts,
2020                  * because an extra constraint step is needed
2021                  */
2022                 do_update_sd2(upd->sd,
2023                               bInitStep, start_th, end_th,
2024                               inputrec->opts.acc, inputrec->opts.nFreeze,
2025                               md->invmass, md->ptype,
2026                               md->cFREEZE, md->cACC, md->cTC,
2027                               state->x, xprime, state->v, force, state->sd_X,
2028                               inputrec->opts.tau_t,
2029                               TRUE, step, inputrec->ld_seed,
2030                               DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2031                 break;
2032             case (eiBD):
2033                 do_update_bd(start_th, end_th, dt,
2034                              inputrec->opts.nFreeze, md->invmass, md->ptype,
2035                              md->cFREEZE, md->cTC,
2036                              state->x, xprime, state->v, force,
2037                              inputrec->bd_fric,
2038                              upd->sd->bd_rf,
2039                              step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->gatindex : NULL);
2040                 break;
2041             case (eiVV):
2042             case (eiVVAK):
2043                 alpha = 1.0 + DIM/((double)inputrec->opts.nrdf[0]); /* assuming barostat coupled to group 0. */
2044                 switch (UpdatePart)
2045                 {
2046                     case etrtVELOCITY1:
2047                     case etrtVELOCITY2:
2048                         do_update_vv_vel(start_th, end_th, dt,
2049                                          inputrec->opts.acc, inputrec->opts.nFreeze,
2050                                          md->invmass, md->ptype,
2051                                          md->cFREEZE, md->cACC,
2052                                          state->v, force,
2053                                          (bNH || bPR), state->veta, alpha);
2054                         break;
2055                     case etrtPOSITION:
2056                         do_update_vv_pos(start_th, end_th, dt,
2057                                          inputrec->opts.nFreeze,
2058                                          md->ptype, md->cFREEZE,
2059                                          state->x, xprime, state->v,
2060                                          (bNH || bPR), state->veta);
2061                         break;
2062                 }
2063                 break;
2064             default:
2065                 gmx_fatal(FARGS, "Don't know how to update coordinates");
2066                 break;
2067         }
2068     }
2069
2070 }
2071
2072
2073 void correct_ekin(FILE *log, int start, int end, rvec v[], rvec vcm, real mass[],
2074                   real tmass, tensor ekin)
2075 {
2076     /*
2077      * This is a debugging routine. It should not be called for production code
2078      *
2079      * The kinetic energy should calculated according to:
2080      *   Ekin = 1/2 m (v-vcm)^2
2081      * However the correction is not always applied, since vcm may not be
2082      * known in time and we compute
2083      *   Ekin' = 1/2 m v^2 instead
2084      * This can be corrected afterwards by computing
2085      *   Ekin = Ekin' + 1/2 m ( -2 v vcm + vcm^2)
2086      * or in hsorthand:
2087      *   Ekin = Ekin' - m v vcm + 1/2 m vcm^2
2088      */
2089     int    i, j, k;
2090     real   m, tm;
2091     rvec   hvcm, mv;
2092     tensor dekin;
2093
2094     /* Local particles */
2095     clear_rvec(mv);
2096
2097     /* Processor dependent part. */
2098     tm = 0;
2099     for (i = start; (i < end); i++)
2100     {
2101         m      = mass[i];
2102         tm    += m;
2103         for (j = 0; (j < DIM); j++)
2104         {
2105             mv[j] += m*v[i][j];
2106         }
2107     }
2108     /* Shortcut */
2109     svmul(1/tmass, vcm, vcm);
2110     svmul(0.5, vcm, hvcm);
2111     clear_mat(dekin);
2112     for (j = 0; (j < DIM); j++)
2113     {
2114         for (k = 0; (k < DIM); k++)
2115         {
2116             dekin[j][k] += vcm[k]*(tm*hvcm[j]-mv[j]);
2117         }
2118     }
2119     pr_rvecs(log, 0, "dekin", dekin, DIM);
2120     pr_rvecs(log, 0, " ekin", ekin, DIM);
2121     fprintf(log, "dekin = %g, ekin = %g  vcm = (%8.4f %8.4f %8.4f)\n",
2122             trace(dekin), trace(ekin), vcm[XX], vcm[YY], vcm[ZZ]);
2123     fprintf(log, "mv = (%8.4f %8.4f %8.4f)\n",
2124             mv[XX], mv[YY], mv[ZZ]);
2125 }
2126
2127 extern gmx_bool update_randomize_velocities(t_inputrec *ir, gmx_int64_t step, const t_commrec *cr,
2128                                             t_mdatoms *md, t_state *state, gmx_update_t upd, gmx_constr_t constr)
2129 {
2130
2131     real rate = (ir->delta_t)/ir->opts.tau_t[0];
2132
2133     if (ir->etc == etcANDERSEN && constr != NULL)
2134     {
2135         gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
2136     }
2137
2138     /* proceed with andersen if 1) it's fixed probability per
2139        particle andersen or 2) it's massive andersen and it's tau_t/dt */
2140     if ((ir->etc == etcANDERSEN) || do_per_step(step, (int)(1.0/rate)))
2141     {
2142         andersen_tcoupl(ir, step, cr, md, state, rate,
2143                         upd->sd->randomize_group, upd->sd->boltzfac);
2144         return TRUE;
2145     }
2146     return FALSE;
2147 }